#include <Digi2Raw2Digi.h>
Public Member Functions | |
virtual void | analyze (const edm::Event &, const edm::EventSetup &) |
virtual void | beginJob () |
template<class Digi > | |
void | compare (const edm::Event &, const edm::EventSetup &) |
Digi2Raw2Digi (const edm::ParameterSet &) | |
virtual void | endJob () |
~Digi2Raw2Digi () | |
Private Attributes | |
DQMStore * | dbe_ |
edm::InputTag | inputTag1_ |
edm::InputTag | inputTag2_ |
MonitorElement * | meStatus |
std::string | outputFile_ |
int | unsuppressed |
Definition at line 14 of file Digi2Raw2Digi.h.
Digi2Raw2Digi::Digi2Raw2Digi | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 220 of file Digi2Raw2Digi.cc.
References DQMStore::book1D(), dbe_, interpolateCardsSimple::histo, meStatus, cppFunctionSkipper::operator, outputFile_, and DQMStore::setCurrentFolder().
: inputTag1_(iConfig.getParameter<edm::InputTag>("digiLabel1")), inputTag2_(iConfig.getParameter<edm::InputTag>("digiLabel2")), outputFile_(iConfig.getUntrackedParameter<std::string>("outputFile")), dbe_(0) { // DQM ROOT output if ( outputFile_.size() != 0 ) { edm::LogInfo("OutputInfo") << " Hcal RecHit Task histograms will be saved to '" << outputFile_.c_str() << "'"; } else { edm::LogInfo("OutputInfo") << " Hcal RecHit Task histograms will NOT be saved"; } // DQM service initialization dbe_ = edm::Service<DQMStore>().operator->(); if ( dbe_ ) { dbe_->setCurrentFolder("Digi2Raw2DigiV/Digi2Raw2DigiTask"); } // const char * sub = hcalselector_.c_str(); char histo[100]; sprintf (histo, "Digi2Raw2Digi_status") ; // bins: 1)full match 2)ID match, not content 3) no match // 4) number of events with diff number of Digis meStatus = dbe_->book1D(histo, histo, 5, 0., 5.); }
Digi2Raw2Digi::~Digi2Raw2Digi | ( | ) |
Definition at line 258 of file Digi2Raw2Digi.cc.
{}
void Digi2Raw2Digi::analyze | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements edm::EDAnalyzer.
Definition at line 261 of file Digi2Raw2Digi.cc.
References iEvent, and unsuppressed.
{ unsuppressed = 0; // std::cout << "=== HBHE ==================" << std::endl; compare<HBHEDataFrame>(iEvent,iSetup); // std::cout << "=== HO ====================" << std::endl; compare<HODataFrame>(iEvent,iSetup); // std::cout << "=== HF ====================" << std::endl; compare<HFDataFrame>(iEvent,iSetup); // std::cout << "=== ZDC ===================" << std::endl; unsuppressed = 1; compare<ZDCDataFrame>(iEvent,iSetup); // std::cout << "=== CASTOR ================" << std::endl; // compare<CastorDataFrame>(iEvent,iSetup); }
void Digi2Raw2Digi::beginJob | ( | void | ) | [virtual] |
void Digi2Raw2Digi::compare | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) |
Definition at line 40 of file Digi2Raw2Digi.cc.
References ecalMGPA::adc(), gather_cfg::cout, HcalDetId::depth(), MonitorElement::Fill(), edm::Event::getByLabel(), i, inputTag1_, inputTag2_, fjr2json::lumi, match(), meStatus, unsuppressed, and HcalZDCDetId::zside().
{ typename edm::Handle<edm::SortedCollection<Digi> > digiCollection1; typename edm::SortedCollection<Digi>::const_iterator digiItr1; typename edm::Handle<edm::SortedCollection<Digi> > digiCollection2; typename edm::SortedCollection<Digi>::const_iterator digiItr2; if(unsuppressed) { // ZDC iEvent.getByLabel ("simHcalUnsuppressedDigis", digiCollection1); } else iEvent.getByLabel (inputTag1_, digiCollection1); iEvent.getByLabel (inputTag2_, digiCollection2); int size1 = 0; int size2 = 0; for (digiItr1=digiCollection1->begin();digiItr1!=digiCollection1->end();digiItr1++) { size1++; } for (digiItr2=digiCollection2->begin();digiItr2!=digiCollection2->end();digiItr2++) { size2++; } //std::cout << "Digi collections size1 = "<< size1 // << " size2 = " << size2 << std::endl; // CYCLE over first DIGI collection ====================================== for (digiItr1=digiCollection1->begin();digiItr1!=digiCollection1->end();digiItr1++) { HcalGenericDetId HcalGenDetId(digiItr1->id()); int tsize = (*digiItr1).size(); int match = 0; if(HcalGenDetId.isHcalZDCDetId()){ //for zdc HcalZDCDetId element(digiItr1->id()); int zside = element.zside(); int section = element.section(); int channel = element.channel(); int gsub = HcalGenDetId.genericSubdet(); if(section==3){// lumi section not reconstructed size2++; match = 1; goto lumi; } //std::cout<< " Zdc genSubdet="<< gsub << " zside=" <<zside // << " section= "<< section << " channel " <<channel // <<std::endl; for (digiItr2=digiCollection2->begin();digiItr2!=digiCollection2->end();digiItr2++) { HcalZDCDetId element2(digiItr2->id()); //int zside2 = element2.zside(); //int section2 = element2.section(); //int channel2 = element2.channel(); //int gsub2 = HcalGenDetId.genericSubdet(); //std::cout<< " Zdc genSubdet="<<gsub2 // << " zside=" <<zside2 // << " section= "<<section2 // <<" channel "<<channel2 // <<std::endl; if(element == element2) { match = 1; int identical = 1; for (int i=0; i<tsize; i++) { double adc = (*digiItr1)[i].adc(); int capid = (*digiItr1)[i].capid(); // std::cout << std::endl << " capid1=" << capid // << " adc1=" << adc // << std::endl; double adc2 = (*digiItr2)[i].adc(); int capid2 = (*digiItr2)[i].capid(); // std::cout << " capid2=" << capid2 << " adc2=" << adc2 // << std::endl; if( capid != capid2 || adc != adc2) { std::cout << "===> PROBLEM !!! gebsubdet=" << gsub << " zside=" <<zside << " section= "<< section << " channel " <<channel << std::endl; std::cout << " capid1["<< i << "]=" << capid << " adc1["<< i << "]=" << adc << " capid2["<< i << "]=" << capid2 << " adc2["<< i << "]=" << adc2 << std::endl; identical = 0; meStatus->Fill(1.); break; } } // end of DataFrames array if(identical) meStatus->Fill(0.); break; // matched HcalZDCID is processed, // go to next (primary collection) cell } } // end of cycle over 2d DIGI collection lumi: if (!match) { meStatus->Fill(2.); std::cout << "===> PROBLEM !!! gsubdet=" << gsub << " zside=" <<zside << " section= "<< section << " channel " <<channel << " HcalZDCId match is not found !!!" << std::endl; } } else{ //for Hcal subdetectors HcalDetId cell(digiItr1->id()); int depth = cell.depth(); int iphi = cell.iphi()-1; int ieta = cell.ieta(); int sub = cell.subdet(); // if(ieta > 0) ieta--; // std::cout << " Cell subdet=" << sub << " ieta=" << ieta // << " inphi=" << iphi << " depth=" << depth << std::endl; // CYCLE over second DIGI collection ====================================== for (digiItr2=digiCollection2->begin();digiItr2!=digiCollection2->end();digiItr2++) { HcalDetId cell2(digiItr2->id()); if( cell == cell2) { match = 1; int identical = 1; for (int i=0; i<tsize; i++) { double adc = (*digiItr1)[i].adc(); int capid = (*digiItr1)[i].capid(); // std::cout << std::endl << " capid1=" << capid // << " adc1=" << adc // << std::endl; double adc2 = (*digiItr2)[i].adc(); int capid2 = (*digiItr2)[i].capid(); // std::cout << " capid2=" << capid2 << " adc2=" << adc2 // << std::endl; if( capid != capid2 || adc != adc2) { std::cout << "===> PROBLEM !!! subdet=" << sub << " ieta=" << ieta << " inphi=" << iphi << " depth=" << depth << std::endl; std::cout << " capid1["<< i << "]=" << capid << " adc1["<< i << "]=" << adc << " capid2["<< i << "]=" << capid2 << " adc2["<< i << "]=" << adc2 << std::endl; identical = 0; meStatus->Fill(1.); break; } } // end of DataFrames array if(identical) meStatus->Fill(0.); break; // matched HcalID is processed, // go to next (primary collection) cell } } // end of cycle over 2d DIGI collection if (!match) { meStatus->Fill(2.); std::cout << "===> PROBLEM !!! subdet=" << sub << " ieta=" << ieta << " inphi=" << iphi << " depth=" << depth << " HcalID match is not found !!!" << std::endl; } } } // end of cycle over 1st DIGI collection if (size1 != size2) { meStatus->Fill(3.); std::cout << "===> PROBLEM !!! Different size of Digi collections : " << size1 << " and " << size2 << std::endl; } }
void Digi2Raw2Digi::endJob | ( | void | ) | [virtual] |
Reimplemented from edm::EDAnalyzer.
Definition at line 254 of file Digi2Raw2Digi.cc.
References dbe_, outputFile_, and DQMStore::save().
{ if ( outputFile_.size() != 0 && dbe_ ) dbe_->save(outputFile_); }
DQMStore* Digi2Raw2Digi::dbe_ [private] |
Definition at line 28 of file Digi2Raw2Digi.h.
Referenced by Digi2Raw2Digi(), and endJob().
edm::InputTag Digi2Raw2Digi::inputTag1_ [private] |
Definition at line 24 of file Digi2Raw2Digi.h.
Referenced by compare().
edm::InputTag Digi2Raw2Digi::inputTag2_ [private] |
Definition at line 25 of file Digi2Raw2Digi.h.
Referenced by compare().
MonitorElement* Digi2Raw2Digi::meStatus [private] |
Definition at line 30 of file Digi2Raw2Digi.h.
Referenced by compare(), and Digi2Raw2Digi().
std::string Digi2Raw2Digi::outputFile_ [private] |
Definition at line 27 of file Digi2Raw2Digi.h.
Referenced by Digi2Raw2Digi(), and endJob().
int Digi2Raw2Digi::unsuppressed [private] |
Definition at line 32 of file Digi2Raw2Digi.h.