00001
00002
00004 #include "HLTriggerOffline/Egamma/interface/EmDQMReco.h"
00005
00007
00009 #include "FWCore/Framework/interface/Frameworkfwd.h"
00010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00011 #include "FWCore/Framework/interface/MakerMacros.h"
00012 #include "DataFormats/HLTReco/interface/TriggerEventWithRefs.h"
00013 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
00014 #include "DataFormats/EgammaCandidates/interface/Electron.h"
00015
00016
00017 #include "DataFormats/L1Trigger/interface/L1EmParticle.h"
00018 #include "DataFormats/L1Trigger/interface/L1EmParticleFwd.h"
00019
00020 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00021 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00022 #include "DataFormats/Common/interface/AssociationMap.h"
00023
00024 #include "DataFormats/Common/interface/Handle.h"
00025 #include "DataFormats/Common/interface/RefToBase.h"
00026 #include "FWCore/ServiceRegistry/interface/Service.h"
00027
00028 #include "FWCore/Utilities/interface/Exception.h"
00029 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
00030 #include "DataFormats/Common/interface/TriggerResults.h"
00031 #include "DataFormats/HLTReco/interface/TriggerObject.h"
00032 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
00034
00036 #include "TFile.h"
00037 #include "TDirectory.h"
00038 #include "TH1F.h"
00039 #include <iostream>
00040 #include <string>
00041 #include <Math/VectorUtil.h>
00042 using namespace ROOT::Math::VectorUtil ;
00043
00044
00046
00048 EmDQMReco::EmDQMReco(const edm::ParameterSet& pset)
00049 {
00050
00051 dbe = edm::Service < DQMStore > ().operator->();
00052 dbe->setVerbose(0);
00053
00054
00056
00058 dirname_="HLT/HLTEgammaValidation/"+pset.getParameter<std::string>("@module_label");
00059 dbe->setCurrentFolder(dirname_);
00060
00061
00062 reqNum = pset.getParameter<unsigned int>("reqNum");
00063 pdgGen = pset.getParameter<int>("pdgGen");
00064 recoEtaAcc = pset.getParameter<double>("genEtaAcc");
00065 recoEtAcc = pset.getParameter<double>("genEtAcc");
00066
00067 plotPtMin = pset.getUntrackedParameter<double>("PtMin",0.);
00068 plotPtMax = pset.getUntrackedParameter<double>("PtMax",1000.);
00069 plotEtaMax = pset.getUntrackedParameter<double>("EtaMax", 2.7);
00070 plotBins = pset.getUntrackedParameter<unsigned int>("Nbins",50);
00071 useHumanReadableHistTitles = pset.getUntrackedParameter<bool>("useHumanReadableHistTitles", false);
00072
00073
00074
00075 recocut_ = pset.getParameter<int>("cutnum");
00076
00077
00078 eventnum = 0;
00079
00080
00081 isHltConfigInitialized_ = false;
00082
00084
00085
00087 std::vector<edm::ParameterSet> filters =
00088 pset.getParameter<std::vector<edm::ParameterSet> >("filters");
00089
00090 int i = 0;
00091 for(std::vector<edm::ParameterSet>::iterator filterconf = filters.begin() ; filterconf != filters.end() ; filterconf++)
00092 {
00093
00094 theHLTCollectionLabels.push_back(filterconf->getParameter<edm::InputTag>("HLTCollectionLabels"));
00095 theHLTOutputTypes.push_back(filterconf->getParameter<int>("theHLTOutputTypes"));
00096
00097 theHLTCollectionHumanNames.push_back(filterconf->getUntrackedParameter<std::string>("HLTCollectionHumanName",theHLTCollectionLabels[i].label()));
00098
00099 std::vector<double> bounds = filterconf->getParameter<std::vector<double> >("PlotBounds");
00100
00101 assert(bounds.size() == 2);
00102 plotBounds.push_back(std::pair<double,double>(bounds[0],bounds[1]));
00103 isoNames.push_back(filterconf->getParameter<std::vector<edm::InputTag> >("IsoCollections"));
00104
00105 assert(isoNames.back().size()>0);
00106 if (isoNames.back().at(0).label()=="none") {
00107 plotiso.push_back(false);
00108 } else {
00109 plotiso.push_back(true);
00110 }
00111 i++;
00112 }
00113
00114
00115 numOfHLTCollectionLabels = theHLTCollectionLabels.size();
00116
00117 }
00118
00119
00120
00124 void EmDQMReco::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup ) {
00125
00126 bool isHltConfigChanged = false;
00127 isHltConfigInitialized_ = hltConfig_.init( iRun, iSetup, "HLT", isHltConfigChanged );
00128
00129 }
00130
00131
00132
00133
00134
00136
00138 void
00139 EmDQMReco::beginJob()
00140 {
00141
00142 dbe->setCurrentFolder(dirname_);
00143
00145
00146
00147
00149
00150 std::string histName="total_eff";
00151 std::string histTitle = "total events passing";
00152
00153
00154 totalreco = dbe->book1D(histName.c_str(),histTitle.c_str(),numOfHLTCollectionLabels+2,0,numOfHLTCollectionLabels+2);
00155 totalreco->setBinLabel(numOfHLTCollectionLabels+1,"Total");
00156 totalreco->setBinLabel(numOfHLTCollectionLabels+2,"Reco");
00157 for (unsigned int u=0; u<numOfHLTCollectionLabels; u++){totalreco->setBinLabel(u+1,theHLTCollectionLabels[u].label().c_str());}
00158
00159 histName="total_eff_RECO_matched";
00160 histTitle="total events passing (Reco matched)";
00161 totalmatchreco = dbe->book1D(histName.c_str(),histTitle.c_str(),numOfHLTCollectionLabels+2,0,numOfHLTCollectionLabels+2);
00162 totalmatchreco->setBinLabel(numOfHLTCollectionLabels+1,"Total");
00163 totalmatchreco->setBinLabel(numOfHLTCollectionLabels+2,"Reco");
00164 for (unsigned int u=0; u<numOfHLTCollectionLabels; u++){totalmatchreco->setBinLabel(u+1,theHLTCollectionLabels[u].label().c_str());}
00165
00166 MonitorElement* tmphisto;
00167 MonitorElement* tmpiso;
00168
00170
00172 std::string pdgIdString;
00173 switch(pdgGen) {
00174 case 11:
00175 pdgIdString="Electron";break;
00176 case 22:
00177 pdgIdString="Photon";break;
00178 default:
00179 pdgIdString="Particle";
00180 }
00181
00182 histName = "reco_et";
00183 histTitle= "E_{T} of " + pdgIdString + "s" ;
00184 etreco = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax);
00185 histName = "reco_eta";
00186 histTitle= "#eta of "+ pdgIdString +"s " ;
00187 etareco = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax);
00188
00189 histName = "reco_et_monpath";
00190 histTitle= "E_{T} of " + pdgIdString + "s monpath" ;
00191 etrecomonpath = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax);
00192 histName = "reco_eta_monpath";
00193 histTitle= "#eta of "+ pdgIdString +"s monpath" ;
00194 etarecomonpath = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax);
00195
00196 histName = "final_et_monpath";
00197 histTitle = "Final Et Monpath";
00198 ethistmonpath = dbe->book1D(histName.c_str(), histTitle.c_str(), plotBins, plotPtMin, plotPtMax);
00199
00200 histName = "final_eta_monpath";
00201 histTitle = "Final Eta Monpath";
00202 etahistmonpath = dbe->book1D(histName.c_str(), histTitle.c_str(), plotBins, -plotEtaMax, plotEtaMax);
00203
00205
00207
00208
00209 std::vector<std::string> HltHistTitle;
00210 if ( theHLTCollectionHumanNames.size() == numOfHLTCollectionLabels && useHumanReadableHistTitles ) {
00211 HltHistTitle = theHLTCollectionHumanNames;
00212 } else {
00213 for (unsigned int i =0; i < numOfHLTCollectionLabels; i++) {
00214 HltHistTitle.push_back(theHLTCollectionLabels[i].label());
00215 }
00216 }
00217
00218 for(unsigned int i = 0; i< numOfHLTCollectionLabels ; i++){
00219
00220 histName = theHLTCollectionLabels[i].label()+"et_all";
00221 histTitle = HltHistTitle[i]+" Et (ALL)";
00222 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax);
00223 ethist.push_back(tmphisto);
00224
00225
00226 histName = theHLTCollectionLabels[i].label()+"eta_all";
00227 histTitle = HltHistTitle[i]+" #eta (ALL)";
00228 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax);
00229 etahist.push_back(tmphisto);
00230
00231
00232 histName = theHLTCollectionLabels[i].label()+"et_RECO_matched";
00233 histTitle = HltHistTitle[i]+" Et (RECO matched)";
00234 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax);
00235 ethistmatchreco.push_back(tmphisto);
00236
00237
00238 histName = theHLTCollectionLabels[i].label()+"eta_RECO_matched";
00239 histTitle = HltHistTitle[i]+" #eta (RECO matched)";
00240 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax);
00241 etahistmatchreco.push_back(tmphisto);
00242
00243
00244
00245 histName = theHLTCollectionLabels[i].label()+"et_RECO_matched_monpath";
00246 histTitle = HltHistTitle[i]+" Et (RECO matched, monpath)";
00247 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax);
00248 ethistmatchrecomonpath.push_back(tmphisto);
00249
00250
00251 histName = theHLTCollectionLabels[i].label()+"eta_RECO_matched_monpath";
00252 histTitle = HltHistTitle[i]+" #eta (RECO matched, monpath)";
00253 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax);
00254 etahistmatchrecomonpath.push_back(tmphisto);
00255
00256
00257 histName = theHLTCollectionLabels[i].label()+"et_reco";
00258 histTitle = HltHistTitle[i]+" Et (reco)";
00259 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax);
00260 histEtOfHltObjMatchToReco.push_back(tmphisto);
00261
00262
00263 histName = theHLTCollectionLabels[i].label()+"eta_reco";
00264 histTitle = HltHistTitle[i]+" eta (reco)";
00265 tmphisto = dbe->book1D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax);
00266 histEtaOfHltObjMatchToReco.push_back(tmphisto);
00267
00268 if (!plotiso[i]) {
00269 tmpiso = NULL;
00270 etahistiso.push_back(tmpiso);
00271 ethistiso.push_back(tmpiso);
00272 etahistisomatchreco.push_back(tmpiso);
00273 ethistisomatchreco.push_back(tmpiso);
00274 histEtaIsoOfHltObjMatchToReco.push_back(tmpiso);
00275 histEtIsoOfHltObjMatchToReco.push_back(tmpiso);
00276 } else {
00277
00278 histName = theHLTCollectionLabels[i].label()+"eta_isolation_all";
00279 histTitle = HltHistTitle[i]+" isolation vs #eta (all)";
00280 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00281 etahistiso.push_back(tmpiso);
00282
00283
00284 histName = theHLTCollectionLabels[i].label()+"et_isolation_all";
00285 histTitle = HltHistTitle[i]+" isolation vs Et (all)";
00286 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00287 ethistiso.push_back(tmpiso);
00288
00289
00290 histName = theHLTCollectionLabels[i].label()+"eta_isolation_RECO_matched";
00291 histTitle = HltHistTitle[i]+" isolation vs #eta (reco matched)";
00292 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00293 etahistisomatchreco.push_back(tmpiso);
00294
00295
00296 histName = theHLTCollectionLabels[i].label()+"et_isolation_RECO_matched";
00297 histTitle = HltHistTitle[i]+" isolation vs Et (reco matched)";
00298 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00299 ethistisomatchreco.push_back(tmpiso);
00300
00301
00302
00303 histName = theHLTCollectionLabels[i].label()+"eta_isolation_reco";
00304 histTitle = HltHistTitle[i]+" isolation vs #eta (reco)";
00305 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,-plotEtaMax,plotEtaMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00306 histEtaIsoOfHltObjMatchToReco.push_back(tmpiso);
00307
00308
00309
00310 histName = theHLTCollectionLabels[i].label()+"et_isolation_reco";
00311 histTitle = HltHistTitle[i]+" isolation vs Et (reco)";
00312 tmpiso = dbe->book2D(histName.c_str(),histTitle.c_str(),plotBins,plotPtMin,plotPtMax,plotBins,plotBounds[i].first,plotBounds[i].second);
00313 histEtIsoOfHltObjMatchToReco.push_back(tmpiso);
00314
00315 }
00316
00317 }
00318 }
00319
00320
00322
00324 EmDQMReco::~EmDQMReco(){
00325 }
00326
00327
00329
00331 void
00332 EmDQMReco::analyze(const edm::Event & event , const edm::EventSetup& setup)
00333 {
00334
00335
00336 if( !isHltConfigInitialized_ ) return;
00337
00338 eventnum++;
00339 bool plotMonpath = false;
00340 bool plotReco = true;
00341
00342 edm::Handle<edm::View<reco::Candidate> > recoObjects;
00343 edm::Handle<std::vector<reco::SuperCluster> > recoObjectsEB;
00344 edm::Handle<std::vector<reco::SuperCluster> > recoObjectsEE;
00345
00346 if (pdgGen == 11) {
00347
00348 event.getByLabel("gsfElectrons", recoObjects);
00349
00350 if (recoObjects->size() < (unsigned int)recocut_) {
00351
00352 return;
00353 }
00354 } else if (pdgGen == 22) {
00355
00356 event.getByLabel("correctedHybridSuperClusters", recoObjectsEB);
00357 event.getByLabel("correctedMulti5x5SuperClustersWithPreshower", recoObjectsEE);
00358
00359 if (recoObjectsEB->size() + recoObjectsEE->size() < (unsigned int)recocut_) {
00360
00361 return;
00362 }
00363 }
00364
00365 edm::Handle<edm::TriggerResults> HLTR;
00366 event.getByLabel(edm::InputTag("TriggerResults","","HLT"), HLTR);
00367
00372
00373
00374
00375
00376
00377
00378
00379
00380
00381
00382
00383 unsigned int triggerIndex;
00384 triggerIndex = hltConfig_.triggerIndex("HLT_MinBias");
00385
00386
00387 bool isFired = false;
00388 if (triggerIndex < HLTR->size()){
00389 isFired = HLTR->accept(triggerIndex);
00390 }
00391
00392
00393
00394 edm::Handle<trigger::TriggerEventWithRefs> triggerObj;
00395 event.getByLabel("hltTriggerSummaryRAW",triggerObj);
00396 if(!triggerObj.isValid()) {
00397 edm::LogWarning("EmDQMReco") << "RAW-type HLT results not found, skipping event";
00398 return;
00399 }
00400
00402
00403
00405 totalreco->Fill(numOfHLTCollectionLabels+0.5);
00406 totalmatchreco->Fill(numOfHLTCollectionLabels+.5);
00407
00408
00409
00410
00411
00412
00413
00414
00415
00416
00417
00418
00419
00420
00421
00422
00423
00425
00426
00428
00429
00430
00431
00433
00435
00436
00437 std::vector<reco::Particle> sortedReco;
00438 if (plotReco == true) {
00439 if (pdgGen == 11) {
00440 for(edm::View<reco::Candidate>::const_iterator recopart = recoObjects->begin(); recopart != recoObjects->end();recopart++){
00441 reco::Particle tmpcand( recopart->charge(), recopart->p4(), recopart->vertex(),recopart->pdgId(),recopart->status() );
00442 sortedReco.push_back(tmpcand);
00443 }
00444 }
00445 else if (pdgGen == 22) {
00446 for(std::vector<reco::SuperCluster>::const_iterator recopart2 = recoObjectsEB->begin(); recopart2 != recoObjectsEB->end();recopart2++){
00447 float en = recopart2->energy();
00448 float er = sqrt(pow(recopart2->x(),2) + pow(recopart2->y(),2) + pow(recopart2->z(),2) );
00449 float px = recopart2->energy()*recopart2->x()/er;
00450 float py = recopart2->energy()*recopart2->y()/er;
00451 float pz = recopart2->energy()*recopart2->z()/er;
00452 reco::Candidate::LorentzVector thisLV(px,py,pz,en);
00453 reco::Particle tmpcand( 0, thisLV, math::XYZPoint(0.,0.,0.), 22, 1 );
00454 sortedReco.push_back(tmpcand);
00455 }
00456 for(std::vector<reco::SuperCluster>::const_iterator recopart2 = recoObjectsEE->begin(); recopart2 != recoObjectsEE->end();recopart2++){
00457 float en = recopart2->energy();
00458 float er = sqrt(pow(recopart2->x(),2) + pow(recopart2->y(),2) + pow(recopart2->z(),2) );
00459 float px = recopart2->energy()*recopart2->x()/er;
00460 float py = recopart2->energy()*recopart2->y()/er;
00461 float pz = recopart2->energy()*recopart2->z()/er;
00462 reco::Candidate::LorentzVector thisLV(px,py,pz,en);
00463 reco::Particle tmpcand( 0, thisLV, math::XYZPoint(0.,0.,0.), 22, 1 );
00464 sortedReco.push_back(tmpcand);
00465 }
00466 }
00467
00468 std::sort(sortedReco.begin(),sortedReco.end(),pTComparator_ );
00469
00470
00471
00472
00473
00474 sortedReco.erase(sortedReco.begin()+recocut_,sortedReco.end());
00475
00476 for (unsigned int i = 0 ; i < recocut_ ; i++ ) {
00477 etreco ->Fill( sortedReco[i].et() );
00478 etareco->Fill( sortedReco[i].eta() );
00479
00480 if (isFired) {
00481 etrecomonpath->Fill( sortedReco[i].et() );
00482 etarecomonpath->Fill( sortedReco[i].eta() );
00483 plotMonpath = true;
00484 }
00485
00486 }
00487
00488 if (recocut_ >= reqNum) totalreco->Fill(numOfHLTCollectionLabels+1.5);
00489 if (recocut_ >= reqNum) totalmatchreco->Fill(numOfHLTCollectionLabels+1.5);
00490
00491 }
00492
00493
00494
00495
00497
00499 for(unsigned int n=0; n < numOfHLTCollectionLabels ; n++) {
00500
00501
00502 switch(theHLTOutputTypes[n])
00503 {
00504 case trigger::TriggerL1NoIsoEG:
00505 fillHistos<l1extra::L1EmParticleCollection>(triggerObj,event,n, sortedReco, plotReco, plotMonpath);break;
00506 case trigger::TriggerL1IsoEG:
00507 fillHistos<l1extra::L1EmParticleCollection>(triggerObj,event,n, sortedReco, plotReco, plotMonpath);break;
00508 case trigger::TriggerPhoton:
00509 fillHistos<reco::RecoEcalCandidateCollection>(triggerObj,event,n, sortedReco, plotReco, plotMonpath);break;
00510 case trigger::TriggerElectron:
00511 fillHistos<reco::ElectronCollection>(triggerObj,event,n, sortedReco, plotReco, plotMonpath);break;
00512 case trigger::TriggerCluster:
00513 fillHistos<reco::RecoEcalCandidateCollection>(triggerObj,event,n, sortedReco, plotReco, plotMonpath);break;
00514 default:
00515 throw(cms::Exception("Release Validation Error") << "HLT output type not implemented: theHLTOutputTypes[n]" );
00516 }
00517 }
00518 }
00519
00520
00522
00523
00525 template <class T> void EmDQMReco::fillHistos(edm::Handle<trigger::TriggerEventWithRefs>& triggerObj,const edm::Event& iEvent ,unsigned int n, std::vector<reco::Particle>& sortedReco, bool plotReco, bool plotMonpath)
00526 {
00527 std::vector<edm::Ref<T> > recoecalcands;
00528 if ( ( triggerObj->filterIndex(theHLTCollectionLabels[n])>=triggerObj->size() )){
00529 return;
00530 }
00531
00533
00535 triggerObj->getObjects(triggerObj->filterIndex(theHLTCollectionLabels[n]),theHLTOutputTypes[n],recoecalcands);
00536
00537
00538 if (theHLTOutputTypes[n] == trigger::TriggerL1NoIsoEG){
00539 std::vector<edm::Ref<T> > isocands;
00540 triggerObj->getObjects(triggerObj->filterIndex(theHLTCollectionLabels[n]),trigger::TriggerL1IsoEG,isocands);
00541 if (isocands.size()>0)
00542 {
00543 for (unsigned int i=0; i < isocands.size(); i++)
00544 recoecalcands.push_back(isocands[i]);
00545 }
00546 }
00547
00548
00549 if (recoecalcands.size() < 1){
00550 return;
00551 }
00552
00553
00554 if (recoecalcands.size() >= reqNum )
00555 totalreco->Fill(n+0.5);
00556
00557
00559
00560
00562 for (unsigned int j=0; j<recoecalcands.size(); j++){
00563 if(!( recoecalcands.at(j).isAvailable())){
00564 edm::LogError("EmDQMReco") << "Event content inconsistent: TriggerEventWithRefs contains invalid Refs" << std::endl << "invalid refs for: " << theHLTCollectionLabels[n].label();
00565 return;
00566 }
00567 }
00568
00570
00571
00573
00574
00575 for (unsigned int i=0; i<recoecalcands.size(); i++) {
00576
00577 ethist[n] ->Fill(recoecalcands[i]->et() );
00578 etahist[n]->Fill(recoecalcands[i]->eta() );
00579
00581
00582
00584 if ( n+1 < numOfHLTCollectionLabels ) {
00585 if (plotiso[n+1]) {
00586 for (unsigned int j = 0 ; j < isoNames[n+1].size() ;j++ ){
00587 edm::Handle<edm::AssociationMap<edm::OneToValue< T , float > > > depMap;
00588 iEvent.getByLabel(isoNames[n+1].at(j),depMap);
00589 if (depMap.isValid()){
00590 typename edm::AssociationMap<edm::OneToValue< T , float > >::const_iterator mapi = depMap->find(recoecalcands[i]);
00591 if (mapi!=depMap->end()){
00592 etahistiso[n+1]->Fill(recoecalcands[i]->eta(),mapi->val);
00593 ethistiso[n+1]->Fill(recoecalcands[i]->et(),mapi->val);
00594 }
00595 }
00596 }
00597 }
00598 }
00599 }
00600
00602
00603
00605 if (plotReco == true) {
00606 for (unsigned int i=0; i < recocut_; i++) {
00607 math::XYZVector currentRecoParticleMomentum = sortedReco[i].momentum();
00608
00609
00610 float closestRecoDeltaR = 1000.;
00611 int closestRecoEcalCandIndex = -1;
00612 for (unsigned int j=0; j<recoecalcands.size(); j++) {
00613 float deltaR = DeltaR(recoecalcands[j]->momentum(),currentRecoParticleMomentum);
00614
00615 if (deltaR < closestRecoDeltaR) {
00616 closestRecoDeltaR = deltaR;
00617 closestRecoEcalCandIndex = j;
00618 }
00619 }
00620
00621
00622
00623 if ( closestRecoEcalCandIndex >= 0 ) {
00624 histEtOfHltObjMatchToReco[n] ->Fill( recoecalcands[closestRecoEcalCandIndex]->et() );
00625 histEtaOfHltObjMatchToReco[n]->Fill( recoecalcands[closestRecoEcalCandIndex]->eta() );
00626
00627
00628 if (n+1 < numOfHLTCollectionLabels){
00629 if (plotiso[n+1] ){
00630 for (unsigned int j = 0 ; j < isoNames[n+1].size() ;j++ ){
00631 edm::Handle<edm::AssociationMap<edm::OneToValue< T , float > > > depMap;
00632 iEvent.getByLabel(isoNames[n+1].at(j),depMap);
00633 if (depMap.isValid()){
00634 typename edm::AssociationMap<edm::OneToValue< T , float > >::const_iterator mapi = depMap->find(recoecalcands[closestRecoEcalCandIndex]);
00635 if (mapi!=depMap->end()) {
00636 histEtaIsoOfHltObjMatchToReco[n+1]->Fill( recoecalcands[closestRecoEcalCandIndex]->eta(),mapi->val);
00637 histEtIsoOfHltObjMatchToReco[n+1] ->Fill( recoecalcands[closestRecoEcalCandIndex]->et(), mapi->val);
00638 }
00639 }
00640 }
00641 }
00642 }
00643 }
00644 }
00645
00647
00649 unsigned int mtachedRecoParts = 0;
00650 float minrecodist=0.3;
00651 if(n==0) minrecodist=0.5;
00652 for(unsigned int i =0; i < recocut_; i++){
00653
00654 bool matchThis= false;
00655 math::XYZVector candDir=sortedReco[i].momentum();
00656 unsigned int closest = 0;
00657 double closestDr = 1000.;
00658 for(unsigned int trigOb = 0 ; trigOb < recoecalcands.size(); trigOb++){
00659 double dr = DeltaR(recoecalcands[trigOb]->momentum(),candDir);
00660 if (dr < closestDr) {
00661 closestDr = dr;
00662 closest = trigOb;
00663 }
00664 if (closestDr > minrecodist) {
00665 closest = -1;
00666 } else {
00667 mtachedRecoParts++;
00668 matchThis = true;
00669 }
00670 }
00671 if ( !matchThis ) continue;
00672
00673 ethistmatchreco[n] ->Fill( sortedReco[i].et() );
00674 etahistmatchreco[n]->Fill( sortedReco[i].eta() );
00675 if (plotMonpath) {
00676 ethistmatchrecomonpath[n]->Fill( sortedReco[i].et() );
00677 etahistmatchrecomonpath[n]->Fill( sortedReco[i].eta() );
00678 }
00680
00681
00683 if (n+1 < numOfHLTCollectionLabels){
00684 if (plotiso[n+1] ){
00685 for (unsigned int j = 0 ; j < isoNames[n+1].size() ;j++ ){
00686 edm::Handle<edm::AssociationMap<edm::OneToValue< T , float > > > depMapReco;
00687 iEvent.getByLabel(isoNames[n+1].at(j),depMapReco);
00688 if (depMapReco.isValid()){
00689 typename edm::AssociationMap<edm::OneToValue< T , float > >::const_iterator mapi = depMapReco->find(recoecalcands[closest]);
00690 if (mapi!=depMapReco->end()){
00691 etahistisomatchreco[n+1]->Fill(sortedReco[i].eta(),mapi->val);
00692 ethistisomatchreco[n+1]->Fill(sortedReco[i].et(),mapi->val);
00693 }
00694 }
00695 }
00696 }
00697 }
00698 }
00699
00700 if (mtachedRecoParts >= reqNum )
00701 totalmatchreco->Fill(n+0.5);
00702 }
00703
00704 }
00705
00706
00708
00710 void EmDQMReco::endJob(){
00711
00712 }
00713
00714 DEFINE_FWK_MODULE(EmDQMReco);