#include <Analyzer_minbias.h>
Public Member Functions | |
virtual void | analyze (const edm::Event &, const edm::EventSetup &) |
Analyzer_minbias (const edm::ParameterSet &) | |
virtual void | beginJob () |
virtual void | beginRun (const edm::Run &r, const edm::EventSetup &iSetup) |
virtual void | endJob () |
virtual void | endRun (const edm::Run &r, const edm::EventSetup &iSetup) |
~Analyzer_minbias () | |
Private Attributes | |
int | depth |
float | eta |
std::string | fOutputFileName |
TH1F * | hbheNoiseE |
edm::InputTag | hbherecoMB |
edm::InputTag | hbherecoNoise |
TH1F * | hbheSignalE |
std::string | hcalfile_ |
TH1F * | hCalo1 [73][43] |
TH1F * | hCalo1mom2 [73][43] |
TH1F * | hCalo2 [73][43] |
TH1F * | hCalo2mom2 [73][43] |
TH1F * | hfNoiseE |
edm::InputTag | hfrecoMB |
edm::InputTag | hfrecoNoise |
TH1F * | hfSignalE |
TH2F * | hHBHEsize_vs_run |
TH2F * | hHFsize_vs_run |
edm::InputTag | horecoMB |
edm::InputTag | horecoNoise |
TFile * | hOutputFile |
int | ieta |
int | iphi |
double | meannoise_min [73][43] |
double | meannoise_pl [73][43] |
float | mom0_Diff |
float | mom0_MB |
float | mom0_Noise |
float | mom1_Diff |
float | mom1_MB |
float | mom1_Noise |
float | mom2_Diff |
float | mom2_MB |
float | mom2_Noise |
float | mom3_Diff |
float | mom3_MB |
float | mom3_Noise |
float | mom4_Diff |
float | mom4_MB |
float | mom4_Noise |
int | mydet |
std::ofstream * | myout_hcal |
int | mysubd |
TTree * | myTree |
double | nevent |
int | nevent_run |
double | noise_min [73][43] |
double | noise_pl [73][43] |
float | occup |
float | phi |
double | theDFFillDetMapMin0 [5][5][73][43] |
double | theDFFillDetMapMin1 [5][5][73][43] |
double | theDFFillDetMapMin2 [5][5][73][43] |
double | theDFFillDetMapPl0 [5][5][73][43] |
double | theDFFillDetMapPl1 [5][5][73][43] |
double | theDFFillDetMapPl2 [5][5][73][43] |
double | theMBFillDetMapMin0 [5][5][73][43] |
double | theMBFillDetMapMin1 [5][5][73][43] |
double | theMBFillDetMapMin2 [5][5][73][43] |
double | theMBFillDetMapMin4 [5][5][73][43] |
double | theMBFillDetMapPl0 [5][5][73][43] |
double | theMBFillDetMapPl1 [5][5][73][43] |
double | theMBFillDetMapPl2 [5][5][73][43] |
double | theMBFillDetMapPl4 [5][5][73][43] |
double | theNSFillDetMapMin0 [5][5][73][43] |
double | theNSFillDetMapMin1 [5][5][73][43] |
double | theNSFillDetMapMin2 [5][5][73][43] |
double | theNSFillDetMapMin4 [5][5][73][43] |
double | theNSFillDetMapPl0 [5][5][73][43] |
double | theNSFillDetMapPl1 [5][5][73][43] |
double | theNSFillDetMapPl2 [5][5][73][43] |
double | theNSFillDetMapPl4 [5][5][73][43] |
bool | theRecalib |
Definition at line 51 of file Analyzer_minbias.h.
cms::Analyzer_minbias::Analyzer_minbias | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 38 of file Analyzer_minbias.cc.
References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), hcalLocalRecoNZS_cff::hbherecoMB, ALCARECOHcalCalMinBias_cff::hbherecoNoise, hcalLocalRecoNZS_cff::hfrecoMB, ALCARECOHcalCalMinBias_cff::hfrecoNoise, hcalLocalRecoNZS_cff::horecoMB, ALCARECOHcalCalMinBias_cff::horecoNoise, i, and j.
{ // get name of output file with histogramms fOutputFileName = iConfig.getUntrackedParameter<string>("HistOutFile"); // get names of modules, producing object collections hbherecoMB = iConfig.getParameter<edm::InputTag>("hbheInputMB"); horecoMB = iConfig.getParameter<edm::InputTag>("hoInputMB"); hfrecoMB = iConfig.getParameter<edm::InputTag>("hfInputMB"); hbherecoNoise = iConfig.getParameter<edm::InputTag>("hbheInputNoise"); horecoNoise = iConfig.getParameter<edm::InputTag>("hoInputNoise"); hfrecoNoise = iConfig.getParameter<edm::InputTag>("hfInputNoise"); theRecalib = iConfig.getParameter<bool>("Recalib"); // // for(int i=0; i<73; i++) { for(int j=0; j<43; j++) { noise_min[i][j] = 0.; noise_pl[i][j] = 0.; } } // // }
cms::Analyzer_minbias::~Analyzer_minbias | ( | ) |
Definition at line 69 of file Analyzer_minbias.cc.
{ // do anything here that needs to be done at desctruction time // (e.g. close files, deallocate resources etc.) }
void cms::Analyzer_minbias::analyze | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 351 of file Analyzer_minbias.cc.
References abs, gather_cfg::cout, HcalDetId::depth(), CaloRecHit::energy(), HcalObjRepresent::Fill(), edm::EventSetup::get(), edm::Event::getByLabel(), L1GtTriggerMenu::gtAlgorithmMap(), hcalLocalRecoNZS_cff::hbherecoMB, ALCARECOHcalCalMinBias_cff::hbherecoNoise, hcalLocalRecoNZS_cff::hfrecoMB, ALCARECOHcalCalMinBias_cff::hfrecoNoise, i, HcalDetId::ieta(), HcalDetId::iphi(), j, gen::k, prof2calltree::l, LogDebug, nevent, funct::pow(), edm::ESHandle< T >::product(), DetId::rawId(), and edm::Event::run().
{ std::cout<<" Start Analyzer_minbias::analyze "<<nevent<<std::endl; nevent++; nevent_run++; using namespace edm; float rnnum = (float)iEvent.run(); /* std::vector<Provenance const*> theProvenance; iEvent.getAllProvenance(theProvenance); for( std::vector<Provenance const*>::const_iterator ip = theProvenance.begin(); ip != theProvenance.end(); ip++) { cout<<" Print all module/label names "<<(**ip).moduleName()<<" "<<(**ip).moduleLabel()<< " "<<(**ip).productInstanceName()<<endl; } */ edm::ESHandle<L1GtTriggerMenu> menuRcd; iSetup.get<L1GtTriggerMenuRcd>().get(menuRcd) ; const L1GtTriggerMenu* menu = menuRcd.product(); const AlgorithmMap& bitMap = menu->gtAlgorithmMap(); edm::Handle<L1GlobalTriggerReadoutRecord> gtRecord; iEvent.getByLabel("gtDigisAlCaMB", gtRecord); if (!gtRecord.isValid()) { // LogDebug("L1GlobalTriggerRecordProducer") // << "\n\n Error: no L1GlobalTriggerReadoutRecord found with input tag " // << m_l1GtReadoutRecord // << "\n Returning empty L1GlobalTriggerRecord.\n\n" // << std::endl; cout<<" No L1 trigger record "<<endl; } else { const DecisionWord dWord = gtRecord->decisionWord(); for (CItAlgo itAlgo = bitMap.begin(); itAlgo != bitMap.end(); itAlgo++) { bool decision=menu->gtAlgorithmResult(itAlgo->first,dWord); if(decision == 1) std::cout<<" Trigger "<<itAlgo->first<<" "<<decision<<std::endl; } } const HcalRespCorrs* myRecalib=0; if( theRecalib ) { // Radek: edm::ESHandle <HcalRespCorrs> recalibCorrs; iSetup.get<HcalRespCorrsRcd>().get("recalibrate",recalibCorrs); myRecalib = recalibCorrs.product(); // end } // theRecalib // Noise part for HB HE double tmpNSFillDetMapPl1[5][5][73][43]; double tmpNSFillDetMapMin1[5][5][73][43]; for (int i=0; i<5;i++) { for (int j=0; j<5;j++) { for (int k=0; k<73;k++) { for (int l=0; l<43;l++) { tmpNSFillDetMapPl1[i][j][k][l] = 0.; tmpNSFillDetMapMin1[i][j][k][l] = 0.; } } } } edm::Handle<HBHERecHitCollection> hbheNormal; iEvent.getByLabel("hbhereco", hbheNormal); if(!hbheNormal.isValid()){ cout<<" hbheNormal failed "<<endl; } else { cout<<" The size of the normal collection "<<hbheNormal->size()<<endl; } edm::Handle<HBHERecHitCollection> hbheNS; iEvent.getByLabel(hbherecoNoise, hbheNS); if(!hbheNS.isValid()){ LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl; cout<<" No HBHE MS "<<endl; return ; } const HBHERecHitCollection HithbheNS = *(hbheNS.product()); cout<<" HBHE NS size of collection "<<HithbheNS.size()<<endl; hHBHEsize_vs_run->Fill(rnnum,(float)HithbheNS.size()); if(HithbheNS.size()!= 5184) { cout<<" HBHE problem "<<rnnum<<" "<<HithbheNS.size()<<endl; return; } edm::Handle<HBHERecHitCollection> hbheMB; iEvent.getByLabel(hbherecoMB, hbheMB); if(!hbheMB.isValid()){ LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl; cout<<" No HBHE MB"<<endl; return ; } const HBHERecHitCollection HithbheMB = *(hbheMB.product()); cout<<" HBHE MB size of collection "<<HithbheMB.size()<<endl; if(HithbheMB.size()!= 5184) { cout<<" HBHE problem "<<rnnum<<" "<<HithbheMB.size()<<endl; return; } edm::Handle<HFRecHitCollection> hfNS; iEvent.getByLabel(hfrecoNoise, hfNS); if(!hfNS.isValid()){ LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl; cout<<" No HF NS "<<endl; return ; } const HFRecHitCollection HithfNS = *(hfNS.product()); cout<<" HFE NS size of collection "<<HithfNS.size()<<endl; hHFsize_vs_run->Fill(rnnum,(float)HithfNS.size()); if(HithfNS.size()!= 1728) { cout<<" HF problem "<<rnnum<<" "<<HithfNS.size()<<endl; return; } edm::Handle<HFRecHitCollection> hfMB; iEvent.getByLabel(hfrecoMB, hfMB); if(!hfMB.isValid()){ LogDebug("") << "HcalCalibAlgos: Error! can't get hbhe product!" << std::endl; cout<<" No HBHE MB"<<endl; return ; } const HFRecHitCollection HithfMB = *(hfMB.product()); cout<<" HF MB size of collection "<<HithfMB.size()<<endl; if(HithfMB.size()!= 1728) { cout<<" HF problem "<<rnnum<<" "<<HithfMB.size()<<endl; return; } for(HBHERecHitCollection::const_iterator hbheItr=HithbheNS.begin(); hbheItr!=HithbheNS.end(); hbheItr++) { // Recalibration of energy float icalconst=1.; DetId mydetid = hbheItr->id().rawId(); if( theRecalib ) icalconst=myRecalib->getValues(mydetid)->getValue(); HBHERecHit aHit(hbheItr->id(),hbheItr->energy()*icalconst,hbheItr->time()); double energyhit = aHit.energy(); DetId id = (*hbheItr).detid(); HcalDetId hid=HcalDetId(id); int mysu = ((hid).rawId()>>25)&0x7; if( hid.ieta() > 0 ) { theNSFillDetMapPl0[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theNSFillDetMapPl0[mysu][hid.depth()][hid.iphi()][hid.ieta()]+ 1.; theNSFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theNSFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()]+energyhit; theNSFillDetMapPl2[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theNSFillDetMapPl2[mysu][hid.depth()][hid.iphi()][hid.ieta()]+pow(energyhit,2); theNSFillDetMapPl4[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theNSFillDetMapPl4[mysu][hid.depth()][hid.iphi()][hid.ieta()]+pow(energyhit,4); tmpNSFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()] = energyhit; } else { theNSFillDetMapMin0[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theNSFillDetMapMin0[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+ 1.; theNSFillDetMapMin1[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theNSFillDetMapMin1[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+energyhit; theNSFillDetMapMin2[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theNSFillDetMapMin2[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+pow(energyhit,2); theNSFillDetMapMin4[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theNSFillDetMapMin4[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+pow(energyhit,4); tmpNSFillDetMapMin1[mysu][hid.depth()][hid.iphi()][hid.ieta()] = energyhit; } if(hid.depth() == 1) { hbheNoiseE->Fill(energyhit); if(energyhit<-2.) std::cout<<" Run "<<rnnum<<" ieta,iphi "<<hid.ieta()<<" "<<hid.iphi()<<energyhit<<std::endl; // if( hid.ieta() > 0 ) { // hCalo1[hid.iphi()][hid.ieta()]->Fill(energyhit-noise_pl[hid.iphi()][hid.ieta()]); // hCalo1mom2[hid.iphi()][hid.ieta()]->Fill(pow(energyhit,2)); // } else { // hCalo2[hid.iphi()][abs(hid.ieta())]->Fill(energyhit-noise_min[hid.iphi()][abs(hid.ieta())]); // hCalo2mom2[hid.iphi()][abs(hid.ieta())]->Fill(pow(energyhit,2)); // } // eta><0 } // depth=1 } // HBHE_NS // Signal part for HB HE for(HBHERecHitCollection::const_iterator hbheItr=HithbheMB.begin(); hbheItr!=HithbheMB.end(); hbheItr++) { // Recalibration of energy float icalconst=1.; DetId mydetid = hbheItr->id().rawId(); if( theRecalib ) icalconst=myRecalib->getValues(mydetid)->getValue(); HBHERecHit aHit(hbheItr->id(),hbheItr->energy()*icalconst,hbheItr->time()); double energyhit = aHit.energy(); DetId id = (*hbheItr).detid(); HcalDetId hid=HcalDetId(id); int mysu = ((hid).rawId()>>25)&0x7; if( hid.ieta() > 0 ) { theMBFillDetMapPl0[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theMBFillDetMapPl0[mysu][hid.depth()][hid.iphi()][hid.ieta()]+ 1.; theMBFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theMBFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()]+energyhit; theMBFillDetMapPl2[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theMBFillDetMapPl2[mysu][hid.depth()][hid.iphi()][hid.ieta()]+pow(energyhit,2); theMBFillDetMapPl4[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theMBFillDetMapPl4[mysu][hid.depth()][hid.iphi()][hid.ieta()]+pow(energyhit,4); float mydiff = energyhit - tmpNSFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()]; theDFFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theDFFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()]+mydiff; theDFFillDetMapPl2[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theDFFillDetMapPl2[mysu][hid.depth()][hid.iphi()][hid.ieta()]+pow(mydiff,2); } else { theMBFillDetMapMin0[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theMBFillDetMapMin0[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+ 1.; theMBFillDetMapMin1[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theMBFillDetMapMin1[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+energyhit; theMBFillDetMapMin2[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theMBFillDetMapMin2[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+pow(energyhit,2); theMBFillDetMapMin4[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theMBFillDetMapMin4[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+pow(energyhit,4); float mydiff = energyhit - tmpNSFillDetMapMin1[mysu][hid.depth()][hid.iphi()][hid.ieta()]; theDFFillDetMapMin1[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theDFFillDetMapMin1[mysu][hid.depth()][hid.iphi()][hid.ieta()]+mydiff; theDFFillDetMapMin2[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theDFFillDetMapMin2[mysu][hid.depth()][hid.iphi()][hid.ieta()]+pow(mydiff,2); } if(hid.depth() == 1) { hbheSignalE->Fill(energyhit); if( hid.ieta() > 0 ) { hCalo1[hid.iphi()][hid.ieta()]->Fill(energyhit); hCalo1mom2[hid.iphi()][hid.ieta()]->Fill(pow(energyhit,2)); } else { hCalo2[hid.iphi()][abs(hid.ieta())]->Fill(energyhit); hCalo2mom2[hid.iphi()][abs(hid.ieta())]->Fill(pow(energyhit,2)); } // eta><0 } // depth=1 } // HBHE_MB // HF for(HFRecHitCollection::const_iterator hbheItr=HithfNS.begin(); hbheItr!=HithfNS.end(); hbheItr++) { // Recalibration of energy float icalconst=1.; DetId mydetid = hbheItr->id().rawId(); if( theRecalib ) icalconst=myRecalib->getValues(mydetid)->getValue(); HFRecHit aHit(hbheItr->id(),hbheItr->energy()*icalconst,hbheItr->time()); double energyhit = aHit.energy(); // // Remove PMT hits // if(fabs(energyhit) > 40. ) continue; DetId id = (*hbheItr).detid(); HcalDetId hid=HcalDetId(id); int mysu = ((hid).rawId()>>25)&0x7; if( hid.ieta() > 0 ) { theNSFillDetMapPl0[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theNSFillDetMapPl0[mysu][hid.depth()][hid.iphi()][hid.ieta()]+ 1.; theNSFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theNSFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()]+energyhit; theNSFillDetMapPl2[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theNSFillDetMapPl2[mysu][hid.depth()][hid.iphi()][hid.ieta()]+pow(energyhit,2); theNSFillDetMapPl4[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theNSFillDetMapPl4[mysu][hid.depth()][hid.iphi()][hid.ieta()]+pow(energyhit,4); tmpNSFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()] = energyhit; } else { theNSFillDetMapMin0[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theNSFillDetMapMin0[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+ 1.; theNSFillDetMapMin1[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theNSFillDetMapMin1[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+energyhit; theNSFillDetMapMin2[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theNSFillDetMapMin2[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+pow(energyhit,2); theNSFillDetMapMin4[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theNSFillDetMapMin4[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+pow(energyhit,4); tmpNSFillDetMapMin1[mysu][hid.depth()][hid.iphi()][hid.ieta()] = energyhit; } if(hid.depth() == 1) { hfNoiseE->Fill(energyhit); //if( hid.ieta() > 0 ) { // hCalo1[hid.iphi()][hid.ieta()]->Fill(energyhit-noise_pl[hid.iphi()][hid.ieta()]); // hCalo1mom2[hid.iphi()][hid.ieta()]->Fill(pow(energyhit,2)); //} else { // hCalo2[hid.iphi()][abs(hid.ieta())]->Fill(energyhit-noise_min[hid.iphi()][abs(hid.ieta())]); // hCalo2mom2[hid.iphi()][abs(hid.ieta())]->Fill(pow(energyhit,2)); //} // eta><0 } // depth=1 } // HBHE_NS // Signal part for HB HE for(HFRecHitCollection::const_iterator hbheItr=HithfMB.begin(); hbheItr!=HithfMB.end(); hbheItr++) { // Recalibration of energy float icalconst=1.; DetId mydetid = hbheItr->id().rawId(); if( theRecalib ) icalconst=myRecalib->getValues(mydetid)->getValue(); HFRecHit aHit(hbheItr->id(),hbheItr->energy()*icalconst,hbheItr->time()); double energyhit = aHit.energy(); // // Remove PMT hits // if(fabs(energyhit) > 40. ) continue; DetId id = (*hbheItr).detid(); HcalDetId hid=HcalDetId(id); int mysu = ((hid).rawId()>>25)&0x7; if( hid.ieta() > 0 ) { theMBFillDetMapPl0[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theMBFillDetMapPl0[mysu][hid.depth()][hid.iphi()][hid.ieta()]+ 1.; theMBFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theMBFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()]+energyhit; theMBFillDetMapPl2[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theMBFillDetMapPl2[mysu][hid.depth()][hid.iphi()][hid.ieta()]+pow(energyhit,2); theMBFillDetMapPl4[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theMBFillDetMapPl4[mysu][hid.depth()][hid.iphi()][hid.ieta()]+pow(energyhit,4); theDFFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theDFFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()]+energyhit-tmpNSFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()]; theDFFillDetMapPl2[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theDFFillDetMapPl2[mysu][hid.depth()][hid.iphi()][hid.ieta()]+pow((energyhit-tmpNSFillDetMapPl1[mysu][hid.depth()][hid.iphi()][hid.ieta()]),2); } else { theMBFillDetMapMin0[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theMBFillDetMapMin0[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+ 1.; theMBFillDetMapMin1[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theMBFillDetMapMin1[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+energyhit; theMBFillDetMapMin2[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theMBFillDetMapMin2[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+pow(energyhit,2); theMBFillDetMapMin4[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())] = theMBFillDetMapMin4[mysu][hid.depth()][hid.iphi()][abs(hid.ieta())]+pow(energyhit,4); theDFFillDetMapMin1[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theDFFillDetMapMin1[mysu][hid.depth()][hid.iphi()][hid.ieta()]+energyhit-tmpNSFillDetMapMin1[mysu][hid.depth()][hid.iphi()][hid.ieta()]; theDFFillDetMapMin2[mysu][hid.depth()][hid.iphi()][hid.ieta()] = theDFFillDetMapMin2[mysu][hid.depth()][hid.iphi()][hid.ieta()]+pow((energyhit-tmpNSFillDetMapMin1[mysu][hid.depth()][hid.iphi()][hid.ieta()]),2); } if(hid.depth() == 1) { hfSignalE->Fill(energyhit); if( hid.ieta() > 0 ) { hCalo1[hid.iphi()][hid.ieta()]->Fill(energyhit); hCalo1mom2[hid.iphi()][hid.ieta()]->Fill(pow(energyhit,2)); } else { hCalo2[hid.iphi()][abs(hid.ieta())]->Fill(energyhit); hCalo2mom2[hid.iphi()][abs(hid.ieta())]->Fill(pow(energyhit,2)); } // eta><0 } // depth=1 } // HF_MB std::cout<<" Event is finished "<<std::endl; }
void cms::Analyzer_minbias::beginJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 86 of file Analyzer_minbias.cc.
References gather_cfg::cout, eta(), i, j, gen::k, prof2calltree::l, nevent, and phi.
{ hOutputFile = new TFile( fOutputFileName.c_str(), "RECREATE" ) ; myTree = new TTree("RecJet","RecJet Tree"); myTree->Branch("mydet", &mydet, "mydet/I"); myTree->Branch("mysubd", &mysubd, "mysubd/I"); myTree->Branch("depth", &depth, "depth/I"); myTree->Branch("ieta", &ieta, "ieta/I"); myTree->Branch("iphi", &iphi, "iphi/I"); myTree->Branch("eta", &eta, "eta/F"); myTree->Branch("phi", &phi, "phi/F"); myTree->Branch("mom0_MB", &mom0_MB, "mom0_MB/F"); myTree->Branch("mom1_MB", &mom1_MB, "mom1_MB/F"); myTree->Branch("mom2_MB", &mom2_MB, "mom2_MB/F"); myTree->Branch("mom4_MB", &mom4_MB, "mom4_MB/F"); myTree->Branch("mom0_Noise", &mom0_Noise, "mom0_Noise/F"); myTree->Branch("mom1_Noise", &mom1_Noise, "mom1_Noise/F"); myTree->Branch("mom2_Noise", &mom2_Noise, "mom2_Noise/F"); myTree->Branch("mom4_Noise", &mom4_Noise, "mom4_Noise/F"); myTree->Branch("mom0_Diff", &mom0_Diff, "mom0_Diff/F"); myTree->Branch("mom1_Diff", &mom1_Diff, "mom1_Diff/F"); myTree->Branch("mom2_Diff", &mom2_Diff, "mom2_Diff/F"); myTree->Branch("occup", &occup, "occup/F"); std::cout<<" Before ordering Histos "<<std::endl; char str0[15]; char str1[15]; char str10[15]; char str11[15]; int k=0; nevent = 0; // Size of collections hHBHEsize_vs_run = new TH2F("hHBHEsize_vs_run","hHBHEsize_vs_run",500,111500.,112000.,6101,-100.5,6000.5); hHFsize_vs_run = new TH2F("hHFsize_vs_run","hHFsize_vs_run",500,111500.,112000.,6101,-100.5,6000.5); for(int i=1;i<73;i++){ for(int j=1;j<43;j++){ meannoise_pl[i][j] = 0.; meannoise_min[i][j] = 0.; // for(int l=1;l<5;l++){ k = i*1000+j; sprintf(str0,"mpl%d",k); sprintf(str1,"mmin%d",k); sprintf(str10,"vpl%d",k); sprintf(str11,"vmin%d",k); // cout<<" "<<i<<" "<<j<<endl; if( j < 30 ) { // first order moment hCalo1[i][j] = new TH1F(str0, "h0", 320, -10., 10.); hCalo2[i][j] = new TH1F(str1, "h1", 320, -10., 10.); // second order moment hCalo1mom2[i][j] = new TH1F(str10, "h10", 320, 0., 20.); hCalo2mom2[i][j] = new TH1F(str11, "h11", 320, 0., 20.); } else { // HF // first order moment // cout<<" "<<i<<" "<<j<<" "<<k<<endl; if(j < 40) { hCalo1[i][j] = new TH1F(str0, "h0", 320, -10., 10.); hCalo2[i][j] = new TH1F(str1, "h1", 320, -10., 10.); // // second order moment hCalo1mom2[i][j] = new TH1F(str10, "h10", 320, 0., 40.); hCalo2mom2[i][j] = new TH1F(str11, "h11", 320, 0., 40.); } else { hCalo1[i][j] = new TH1F(str0,"h0" , 320, -10., 10.); hCalo2[i][j] = new TH1F(str1, "h1", 320, -10., 10.); // second order moment hCalo1mom2[i][j] = new TH1F(str10, "h10", 320, 0., 120.); hCalo2mom2[i][j] = new TH1F(str11, "h11", 320, 0., 120.); } } // HE/HF boundary // } // l } // j } // i hbheNoiseE = new TH1F("hbheNoiseE","hbheNoiseE", 320, -10., 10.); hfNoiseE = new TH1F("hfNoiseE","hfNoiseE", 320, -10., 10.); hbheSignalE = new TH1F("hbheSignalE","hbheSignalE", 320, -10., 10.); hfSignalE = new TH1F("hfSignalE","hfSignalE", 320, -10., 10.); std::cout<<" After ordering Histos "<<std::endl; std::string ccc = "noise_0.dat"; myout_hcal = new ofstream(ccc.c_str()); if(!myout_hcal) cout << " Output file not open!!! "<<endl; // for (int i=0; i<5;i++) { for (int j=0; j<5;j++) { for (int k=0; k<73;k++) { for (int l=0; l<43;l++) { theMBFillDetMapPl0[i][j][k][l] = 0.; theMBFillDetMapPl1[i][j][k][l] = 0.; theMBFillDetMapPl2[i][j][k][l] = 0.; theMBFillDetMapPl4[i][j][k][l] = 0.; theMBFillDetMapMin0[i][j][k][l] = 0.; theMBFillDetMapMin1[i][j][k][l] = 0.; theMBFillDetMapMin2[i][j][k][l] = 0.; theMBFillDetMapMin4[i][j][k][l] = 0.; theNSFillDetMapPl0[i][j][k][l] = 0.; theNSFillDetMapPl1[i][j][k][l] = 0.; theNSFillDetMapPl2[i][j][k][l] = 0.; theNSFillDetMapPl4[i][j][k][l] = 0.; theNSFillDetMapMin0[i][j][k][l] = 0.; theNSFillDetMapMin1[i][j][k][l] = 0.; theNSFillDetMapMin2[i][j][k][l] = 0.; theNSFillDetMapMin4[i][j][k][l] = 0.; theDFFillDetMapPl0[i][j][k][l] = 0.; theDFFillDetMapPl1[i][j][k][l] = 0.; theDFFillDetMapPl2[i][j][k][l] = 0.; theDFFillDetMapMin0[i][j][k][l] = 0.; theDFFillDetMapMin1[i][j][k][l] = 0.; theDFFillDetMapMin2[i][j][k][l] = 0.; } } } } return ; }
void cms::Analyzer_minbias::beginRun | ( | const edm::Run & | r, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 77 of file Analyzer_minbias.cc.
{ nevent_run = 0; }
void cms::Analyzer_minbias::endJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 244 of file Analyzer_minbias.cc.
References gather_cfg::cout, i, j, gen::k, and prof2calltree::l.
{ int ii=0; for (int i=1; i<5;i++) { for (int j=1; j<5;j++) { for (int k=1; k<73;k++) { for (int l=1; l<43;l++) { if(theMBFillDetMapPl0[i][j][k][l] > 0) { mom0_MB = theMBFillDetMapPl0[i][j][k][l]; mom1_MB = theMBFillDetMapPl1[i][j][k][l]; mom2_MB = theMBFillDetMapPl2[i][j][k][l]; mom4_MB = theMBFillDetMapPl4[i][j][k][l]; mom0_Noise = theNSFillDetMapPl0[i][j][k][l]; mom1_Noise = theNSFillDetMapPl1[i][j][k][l]; mom2_Noise = theNSFillDetMapPl2[i][j][k][l]; mom4_Noise = theNSFillDetMapPl4[i][j][k][l]; mom0_Diff = theDFFillDetMapPl0[i][j][k][l]; mom1_Diff = theDFFillDetMapPl1[i][j][k][l]; mom2_Diff = theDFFillDetMapPl2[i][j][k][l]; mysubd = i; depth = j; ieta = l; iphi = k; cout<<" Result Plus= "<<mysubd<<" "<<ieta<<" "<<iphi<<" mom0 "<<mom0_MB<<" mom1 "<<mom1_MB<<" mom2 "<<mom2_MB<<endl; myTree->Fill(); ii++; } // Pl > 0 if(theMBFillDetMapMin0[i][j][k][l] > 0) { mom0_MB = theMBFillDetMapMin0[i][j][k][l]; mom1_MB = theMBFillDetMapMin1[i][j][k][l]; mom2_MB = theMBFillDetMapMin2[i][j][k][l]; mom4_MB = theMBFillDetMapMin4[i][j][k][l]; mom0_Noise = theNSFillDetMapMin0[i][j][k][l]; mom1_Noise = theNSFillDetMapMin1[i][j][k][l]; mom2_Noise = theNSFillDetMapMin2[i][j][k][l]; mom4_Noise = theNSFillDetMapMin4[i][j][k][l]; mom0_Diff = theDFFillDetMapMin0[i][j][k][l]; mom1_Diff = theDFFillDetMapMin1[i][j][k][l]; mom2_Diff = theDFFillDetMapMin2[i][j][k][l]; mysubd = i; depth = j; ieta = -1*l; iphi = k; cout<<" Result Minus= "<<mysubd<<" "<<ieta<<" "<<iphi<<" mom0 "<<mom0_MB<<" mom1 "<<mom1_MB<<" mom2 "<<mom2_MB<<endl; myTree->Fill(); ii++; } // Min>0 } // ieta } // iphi } // depth } //subd cout<<" Number of cells "<<ii<<endl; hOutputFile->Write(); hOutputFile->cd(); myTree->Write(); hHBHEsize_vs_run->Write() ; hHFsize_vs_run->Write() ; for(int i=1;i<73;i++){ for(int j=1;j<43;j++){ hCalo1[i][j]->Write(); hCalo2[i][j]->Write(); hCalo1mom2[i][j]->Write(); hCalo2mom2[i][j]->Write(); } } hbheNoiseE->Write() ; hfNoiseE->Write() ; hbheSignalE->Write() ; hfSignalE->Write() ; hOutputFile->Close() ; cout<<" File is closed "<<endl; return ; }
void cms::Analyzer_minbias::endRun | ( | const edm::Run & | r, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 81 of file Analyzer_minbias.cc.
References gather_cfg::cout, and edm::RunBase::run().
{ std::cout<<" Runnumber "<<r.run()<<" Nevents "<<nevent_run<<std::endl; }
int cms::Analyzer_minbias::depth [private] |
Definition at line 85 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::eta [private] |
Definition at line 86 of file Analyzer_minbias.h.
std::string cms::Analyzer_minbias::fOutputFileName [private] |
Definition at line 64 of file Analyzer_minbias.h.
TH1F* cms::Analyzer_minbias::hbheNoiseE [private] |
Definition at line 76 of file Analyzer_minbias.h.
Definition at line 126 of file Analyzer_minbias.h.
Definition at line 130 of file Analyzer_minbias.h.
TH1F* cms::Analyzer_minbias::hbheSignalE [private] |
Definition at line 77 of file Analyzer_minbias.h.
std::string cms::Analyzer_minbias::hcalfile_ [private] |
Definition at line 65 of file Analyzer_minbias.h.
TH1F* cms::Analyzer_minbias::hCalo1[73][43] [private] |
Definition at line 72 of file Analyzer_minbias.h.
TH1F* cms::Analyzer_minbias::hCalo1mom2[73][43] [private] |
Definition at line 74 of file Analyzer_minbias.h.
TH1F* cms::Analyzer_minbias::hCalo2[73][43] [private] |
Definition at line 73 of file Analyzer_minbias.h.
TH1F* cms::Analyzer_minbias::hCalo2mom2[73][43] [private] |
Definition at line 75 of file Analyzer_minbias.h.
TH1F* cms::Analyzer_minbias::hfNoiseE [private] |
Definition at line 78 of file Analyzer_minbias.h.
edm::InputTag cms::Analyzer_minbias::hfrecoMB [private] |
Definition at line 128 of file Analyzer_minbias.h.
Definition at line 132 of file Analyzer_minbias.h.
TH1F* cms::Analyzer_minbias::hfSignalE [private] |
Definition at line 79 of file Analyzer_minbias.h.
TH2F* cms::Analyzer_minbias::hHBHEsize_vs_run [private] |
Definition at line 81 of file Analyzer_minbias.h.
TH2F* cms::Analyzer_minbias::hHFsize_vs_run [private] |
Definition at line 82 of file Analyzer_minbias.h.
edm::InputTag cms::Analyzer_minbias::horecoMB [private] |
Definition at line 127 of file Analyzer_minbias.h.
Definition at line 131 of file Analyzer_minbias.h.
TFile* cms::Analyzer_minbias::hOutputFile [private] |
Definition at line 70 of file Analyzer_minbias.h.
int cms::Analyzer_minbias::ieta [private] |
Definition at line 85 of file Analyzer_minbias.h.
int cms::Analyzer_minbias::iphi [private] |
Definition at line 85 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::meannoise_min[73][43] [private] |
Definition at line 93 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::meannoise_pl[73][43] [private] |
Definition at line 93 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom0_Diff [private] |
Definition at line 89 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom0_MB [private] |
Definition at line 87 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom0_Noise [private] |
Definition at line 88 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom1_Diff [private] |
Definition at line 89 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom1_MB [private] |
Definition at line 87 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom1_Noise [private] |
Definition at line 88 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom2_Diff [private] |
Definition at line 89 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom2_MB [private] |
Definition at line 87 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom2_Noise [private] |
Definition at line 88 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom3_Diff [private] |
Definition at line 89 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom3_MB [private] |
Definition at line 87 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom3_Noise [private] |
Definition at line 88 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom4_Diff [private] |
Definition at line 89 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom4_MB [private] |
Definition at line 87 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::mom4_Noise [private] |
Definition at line 88 of file Analyzer_minbias.h.
int cms::Analyzer_minbias::mydet [private] |
Definition at line 85 of file Analyzer_minbias.h.
std::ofstream* cms::Analyzer_minbias::myout_hcal [private] |
Definition at line 66 of file Analyzer_minbias.h.
int cms::Analyzer_minbias::mysubd [private] |
Definition at line 85 of file Analyzer_minbias.h.
TTree* cms::Analyzer_minbias::myTree [private] |
Definition at line 71 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::nevent [private] |
Definition at line 98 of file Analyzer_minbias.h.
int cms::Analyzer_minbias::nevent_run [private] |
Definition at line 84 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::noise_min[73][43] [private] |
Definition at line 94 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::noise_pl[73][43] [private] |
Definition at line 94 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::occup [private] |
Definition at line 87 of file Analyzer_minbias.h.
float cms::Analyzer_minbias::phi [private] |
Definition at line 86 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theDFFillDetMapMin0[5][5][73][43] [private] |
Definition at line 122 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theDFFillDetMapMin1[5][5][73][43] [private] |
Definition at line 123 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theDFFillDetMapMin2[5][5][73][43] [private] |
Definition at line 124 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theDFFillDetMapPl0[5][5][73][43] [private] |
Definition at line 119 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theDFFillDetMapPl1[5][5][73][43] [private] |
Definition at line 120 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theDFFillDetMapPl2[5][5][73][43] [private] |
Definition at line 121 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theMBFillDetMapMin0[5][5][73][43] [private] |
Definition at line 104 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theMBFillDetMapMin1[5][5][73][43] [private] |
Definition at line 105 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theMBFillDetMapMin2[5][5][73][43] [private] |
Definition at line 106 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theMBFillDetMapMin4[5][5][73][43] [private] |
Definition at line 107 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theMBFillDetMapPl0[5][5][73][43] [private] |
Definition at line 99 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theMBFillDetMapPl1[5][5][73][43] [private] |
Definition at line 100 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theMBFillDetMapPl2[5][5][73][43] [private] |
Definition at line 101 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theMBFillDetMapPl4[5][5][73][43] [private] |
Definition at line 102 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theNSFillDetMapMin0[5][5][73][43] [private] |
Definition at line 114 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theNSFillDetMapMin1[5][5][73][43] [private] |
Definition at line 115 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theNSFillDetMapMin2[5][5][73][43] [private] |
Definition at line 116 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theNSFillDetMapMin4[5][5][73][43] [private] |
Definition at line 117 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theNSFillDetMapPl0[5][5][73][43] [private] |
Definition at line 109 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theNSFillDetMapPl1[5][5][73][43] [private] |
Definition at line 110 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theNSFillDetMapPl2[5][5][73][43] [private] |
Definition at line 111 of file Analyzer_minbias.h.
double cms::Analyzer_minbias::theNSFillDetMapPl4[5][5][73][43] [private] |
Definition at line 112 of file Analyzer_minbias.h.
bool cms::Analyzer_minbias::theRecalib [private] |
Definition at line 135 of file Analyzer_minbias.h.