CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/Validation/RecoB/plugins/BTagPerformanceAnalyzerMC.cc

Go to the documentation of this file.
00001 #include "FWCore/Framework/interface/MakerMacros.h"
00002 #include "FWCore/Framework/interface/Frameworkfwd.h"
00003 #include "FWCore/Utilities/interface/Exception.h"
00004 #include "FWCore/Utilities/interface/CodedException.h"
00005 #include "Validation/RecoB/plugins/BTagPerformanceAnalyzerMC.h"
00006 #include "DQMOffline/RecoB/interface/JetTagPlotter.h"
00007 #include "DQMOffline/RecoB/interface/TagInfoPlotterFactory.h"
00008 #include "DataFormats/Common/interface/View.h"
00009 #include "DQMServices/Core/interface/DQMStore.h"
00010 #include "FWCore/ServiceRegistry/interface/Service.h"
00011 
00012 using namespace reco;
00013 using namespace edm;
00014 using namespace std;
00015 using namespace RecoBTag;
00016 //using namespace BTagMCTools;
00017 
00018 typedef std::pair<Jet, reco::JetFlavour> JetWithFlavour;
00019 
00020 BTagPerformanceAnalyzerMC::BTagPerformanceAnalyzerMC(const edm::ParameterSet& pSet) :
00021   partonKinematics(pSet.getParameter< bool >("partonKinematics")),
00022   ptPartonMin(pSet.getParameter<double>("ptPartonMin")),
00023   ptPartonMax(pSet.getParameter<double>("ptPartonMax")),
00024   jetSelector(
00025     pSet.getParameter<double>("etaMin"),
00026     pSet.getParameter<double>("etaMax"),
00027     pSet.getParameter<double>("ptRecJetMin"),
00028     pSet.getParameter<double>("ptRecJetMax"),
00029     0.0, 99999.0,
00030     pSet.getParameter<double>("ratioMin"),
00031     pSet.getParameter<double>("ratioMax")
00032   ),
00033   etaRanges(pSet.getParameter< vector<double> >("etaRanges")),
00034   ptRanges(pSet.getParameter< vector<double> >("ptRanges")),
00035   produceEps(pSet.getParameter< bool >("produceEps")),
00036   producePs(pSet.getParameter< bool >("producePs")),
00037   inputFile(pSet.getParameter<std::string>( "inputfile" )),
00038   update(pSet.getParameter<bool>( "update" )),
00039   allHisto(pSet.getParameter<bool>( "allHistograms" )),
00040   finalize(pSet.getParameter< bool >("finalizePlots")),
00041   finalizeOnly(pSet.getParameter< bool >("finalizeOnly")),
00042   jetMCSrc(pSet.getParameter<edm::InputTag>("jetMCSrc")),
00043   slInfoTag(pSet.getParameter<edm::InputTag>("softLeptonInfo")),
00044   moduleConfig(pSet.getParameter< vector<edm::ParameterSet> >("tagConfig")),
00045   mcPlots_(pSet.getParameter< bool >("mcPlots")),
00046   makeDiffPlots_(pSet.getParameter< bool >("differentialPlots")),
00047   jetCorrector(pSet.getParameter<std::string>("jetCorrection")),
00048   jetMatcher(pSet.getParameter<edm::ParameterSet>("recJetMatching"))
00049 {
00050   double ptRecJetMin = pSet.getParameter<double>("ptRecJetMin");
00051   jetMatcher.setThreshold(0.25 * ptRecJetMin);
00052   switch(pSet.getParameter<unsigned int>("leptonPlots")) {
00053     case 11: electronPlots = true; muonPlots = false; tauPlots = false; break;
00054     case 13: muonPlots = true; electronPlots = false; tauPlots = false; break;
00055     case 15: tauPlots = true; electronPlots = false; tauPlots = false; break;
00056     default: electronPlots = false; muonPlots = false; tauPlots = false;
00057   }
00058   bookHistos(pSet);
00059 }
00060 
00061 void BTagPerformanceAnalyzerMC::bookHistos(const edm::ParameterSet& pSet)
00062 {
00063   //
00064   // Book all histograms.
00065   //
00066 
00067   //if (update) {
00068     //
00069     // append the DQM file ... we should consider this experimental
00070     //    edm::Service<DQMStore>().operator->()->open(std::string((const char *)(inputFile)),"/");
00071     // removed; DQM framework will take care
00072   //}
00073 
00074   // parton p
00075 //   double pPartonMin = 0.0    ;
00076 //   double pPartonMax = 99999.9 ;
00077 
00078 
00079   // iterate over ranges:
00080   const int iEtaStart = -1                   ;  // this will be the inactive one
00081   const int iEtaEnd   = etaRanges.size() - 1 ;
00082   const int iPtStart  = -1                   ;  // this will be the inactive one
00083   const int iPtEnd    = ptRanges.size() - 1  ;
00084   setTDRStyle();
00085 
00086   TagInfoPlotterFactory theFactory;
00087 
00088   for (vector<edm::ParameterSet>::const_iterator iModule = moduleConfig.begin();
00089        iModule != moduleConfig.end(); ++iModule) {
00090 
00091     const string& dataFormatType = iModule->exists("type") ?
00092                                    iModule->getParameter<string>("type") :
00093                                    "JetTag";
00094     if (dataFormatType == "JetTag") {
00095       const InputTag& moduleLabel = iModule->getParameter<InputTag>("label");
00096       const string& folderName    = iModule->getParameter<string>("folder");
00097 
00098       jetTagInputTags.push_back(moduleLabel);
00099       binJetTagPlotters.push_back(vector<JetTagPlotter*>()) ;
00100       // Contains plots for each bin of rapidity and pt.
00101         vector<BTagDifferentialPlot*> * differentialPlotsConstantEta = new vector<BTagDifferentialPlot*> () ;
00102         vector<BTagDifferentialPlot*> * differentialPlotsConstantPt  = new vector<BTagDifferentialPlot*> () ;
00103       if (finalize && mcPlots_ && makeDiffPlots_){
00104         differentialPlots.push_back(vector<BTagDifferentialPlot*>());
00105 
00106         // the constant b-efficiency for the differential plots versus pt and eta
00107         const double& effBConst =
00108                                 iModule->getParameter<edm::ParameterSet>("parameters").getParameter<double>("effBConst");
00109 
00110         // the objects for the differential plots vs. eta,pt for
00111         for ( int iEta = iEtaStart ; iEta < iEtaEnd ; iEta++ ) {
00112           BTagDifferentialPlot * etaConstDifferentialPlot = new BTagDifferentialPlot
00113             (effBConst, BTagDifferentialPlot::constETA, folderName);
00114           differentialPlotsConstantEta->push_back ( etaConstDifferentialPlot );
00115         }
00116 
00117         for ( int iPt = iPtStart ; iPt < iPtEnd ; iPt++ ) {
00118           // differentialPlots for this pt bin
00119           BTagDifferentialPlot * ptConstDifferentialPlot = new BTagDifferentialPlot
00120             (effBConst, BTagDifferentialPlot::constPT, folderName);
00121           differentialPlotsConstantPt->push_back ( ptConstDifferentialPlot );
00122         }
00123       }
00124       // eta loop
00125       for ( int iEta = iEtaStart ; iEta < iEtaEnd ; iEta++ ) {
00126         // pt loop
00127         for ( int iPt = iPtStart ; iPt < iPtEnd ; iPt++ ) {
00128 
00129           const EtaPtBin& etaPtBin = getEtaPtBin(iEta, iPt);
00130 
00131           // Instantiate the genertic b tag plotter
00132           JetTagPlotter *jetTagPlotter = new JetTagPlotter(folderName, etaPtBin,
00133                                                            iModule->getParameter<edm::ParameterSet>("parameters"),mcPlots_,update,finalize);
00134           binJetTagPlotters.back().push_back ( jetTagPlotter ) ;
00135 
00136           // Add to the corresponding differential plotters
00137           if (finalize && mcPlots_ && makeDiffPlots_){
00138             (*differentialPlotsConstantEta)[iEta+1]->addBinPlotter ( jetTagPlotter ) ;
00139             (*differentialPlotsConstantPt )[iPt+1] ->addBinPlotter ( jetTagPlotter ) ;
00140           }
00141         }
00142       }
00143       // the objects for the differential plots vs. eta, pt: collect all from constant eta and constant pt
00144       if (finalize && mcPlots_ && makeDiffPlots_){
00145         differentialPlots.back().reserve(differentialPlotsConstantEta->size()+differentialPlotsConstantPt->size()) ;
00146         differentialPlots.back().insert(differentialPlots.back().end(), differentialPlotsConstantEta->begin(), differentialPlotsConstantEta->end());
00147         differentialPlots.back().insert(differentialPlots.back().end(), differentialPlotsConstantPt->begin(), differentialPlotsConstantPt->end());
00148 
00149         edm::LogInfo("Info")
00150           << "====>>>> ## sizeof differentialPlots = " << differentialPlots.size();
00151 
00152         // the intermediate ones are no longer needed
00153         delete differentialPlotsConstantEta ;
00154         delete differentialPlotsConstantPt  ;
00155       }
00156     } else if(dataFormatType == "TagCorrelation") {
00157         const InputTag& label1 = iModule->getParameter<InputTag>("label1");
00158         const InputTag& label2 = iModule->getParameter<InputTag>("label2");
00159         tagCorrelationInputTags.push_back(std::pair<edm::InputTag, edm::InputTag>(label1, label2));
00160         binTagCorrelationPlotters.push_back(vector<TagCorrelationPlotter*>());
00161 
00162         // eta loop
00163         for ( int iEta = iEtaStart ; iEta != iEtaEnd ; ++iEta) {
00164           // pt loop
00165           for( int iPt = iPtStart ; iPt != iPtEnd ; ++iPt) {
00166             const EtaPtBin& etaPtBin = getEtaPtBin(iEta, iPt);
00167             // Instantiate the generic b tag correlation plotter
00168             TagCorrelationPlotter* tagCorrelationPlotter = new TagCorrelationPlotter(label1.label(), label2.label(), etaPtBin,
00169                                                                                      iModule->getParameter<edm::ParameterSet>("parameters"),
00170                                                                                      mcPlots_, update);
00171             binTagCorrelationPlotters.back().push_back(tagCorrelationPlotter);
00172           }
00173         }
00174     } else {
00175       // tag info retrievel is deferred (needs availability of EventSetup)
00176       const InputTag& moduleLabel = iModule->getParameter<InputTag>("label");
00177       const string& folderName    = iModule->getParameter<string>("folder");
00178 
00179       tagInfoInputTags.push_back(vector<edm::InputTag>());
00180       tiDataFormatType.push_back(dataFormatType);
00181       binTagInfoPlotters.push_back(vector<BaseTagInfoPlotter*>()) ;
00182       // eta loop
00183       for ( int iEta = iEtaStart ; iEta < iEtaEnd ; iEta++ ) {
00184         // pt loop
00185         for ( int iPt = iPtStart ; iPt < iPtEnd ; iPt++ ) {
00186           const EtaPtBin& etaPtBin = getEtaPtBin(iEta, iPt);
00187 
00188           // Instantiate the tagInfo plotter
00189 
00190           BaseTagInfoPlotter *jetTagPlotter = theFactory.buildPlotter(dataFormatType, moduleLabel.label(), 
00191                                      etaPtBin, iModule->getParameter<edm::ParameterSet>("parameters"), folderName, 
00192                                      update, mcPlots_,finalize);
00193           binTagInfoPlotters.back().push_back ( jetTagPlotter ) ;
00194           binTagInfoPlottersToModuleConfig.insert(make_pair(jetTagPlotter, iModule - moduleConfig.begin()));
00195         }
00196       }
00197 
00198       edm::LogInfo("Info")
00199         << "====>>>> ## sizeof differentialPlots = " << differentialPlots.size();
00200     }
00201   }
00202 }
00203 
00204 EtaPtBin BTagPerformanceAnalyzerMC::getEtaPtBin(const int& iEta, const int& iPt)
00205 {
00206   // DEFINE BTagBin:
00207   bool    etaActive_ , ptActive_;
00208   double  etaMin_, etaMax_, ptMin_, ptMax_ ;
00209 
00210   if ( iEta != -1 ) {
00211     etaActive_ = true ;
00212     etaMin_    = etaRanges[iEta]   ;
00213     etaMax_    = etaRanges[iEta+1] ;
00214   }
00215   else {
00216     etaActive_ = false ;
00217     etaMin_    = etaRanges[0]   ;
00218     etaMax_    = etaRanges[etaRanges.size() - 1]   ;
00219   }
00220 
00221   if ( iPt != -1 ) {
00222     ptActive_ = true ;
00223     ptMin_    = ptRanges[iPt]   ;
00224     ptMax_    = ptRanges[iPt+1] ;
00225   }
00226   else {
00227     ptActive_ = false ;
00228     ptMin_    = ptRanges[0]     ;
00229     ptMax_    = ptRanges[ptRanges.size() - 1]   ;
00230   }
00231   return EtaPtBin(etaActive_ , etaMin_ , etaMax_ ,
00232                         ptActive_  , ptMin_  , ptMax_ );
00233 }
00234 
00235 BTagPerformanceAnalyzerMC::~BTagPerformanceAnalyzerMC()
00236 {
00237   for (vector<vector<JetTagPlotter*> >::iterator iJetLabel = binJetTagPlotters.begin();
00238        iJetLabel != binJetTagPlotters.end(); ++iJetLabel) 
00239     for (vector<JetTagPlotter*>::iterator iPlotter = iJetLabel->begin(); iPlotter != iJetLabel->end(); ++iPlotter) 
00240       delete *iPlotter;
00241 
00242   if (finalize  && mcPlots_ && makeDiffPlots_) {
00243     for(vector<vector<BTagDifferentialPlot*> >::iterator iJetLabel = differentialPlots.begin();
00244         iJetLabel != differentialPlots.end(); ++iJetLabel)
00245       for (vector<BTagDifferentialPlot *>::iterator iPlotter = iJetLabel->begin();
00246            iPlotter != iJetLabel->end(); ++ iPlotter) 
00247         delete *iPlotter;
00248   }
00249 
00250   for (vector<vector<TagCorrelationPlotter*> >::iterator iJetLabel = binTagCorrelationPlotters.begin(); 
00251        iJetLabel != binTagCorrelationPlotters.end(); ++iJetLabel) 
00252     for (vector<TagCorrelationPlotter* >::iterator iPlotter = iJetLabel->begin(); iPlotter != iJetLabel->end(); ++iPlotter) 
00253       delete *iPlotter;
00254     
00255   
00256   for (vector<vector<BaseTagInfoPlotter*> >::iterator iJetLabel = binTagInfoPlotters.begin(); 
00257        iJetLabel != binTagInfoPlotters.end(); ++iJetLabel) 
00258     for (vector<BaseTagInfoPlotter*>::iterator iPlotter = iJetLabel->begin(); iPlotter != iJetLabel->end(); ++iPlotter) 
00259       delete *iPlotter;
00260     
00261 }
00262 
00263 void BTagPerformanceAnalyzerMC::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00264 {
00265   eventInitialized = false;
00266 
00267   if (finalizeOnly) return;
00268 
00269   edm::Handle<JetFlavourMatchingCollection> jetMC;
00270   FlavourMap flavours;
00271   LeptonMap leptons;
00272 
00273   iEvent.getByLabel(jetMCSrc, jetMC);
00274   for (JetFlavourMatchingCollection::const_iterator iter = jetMC->begin();
00275        iter != jetMC->end(); ++iter) {
00276     unsigned int fl = std::abs(iter->second.getFlavour());
00277     flavours.insert(std::make_pair(iter->first, fl));
00278     const reco::JetFlavour::Leptons &lep = iter->second.getLeptons();
00279     leptons.insert(std::make_pair(iter->first, lep));
00280   }
00281 
00282   edm::Handle<reco::SoftLeptonTagInfoCollection> infoHandle;
00283   iEvent.getByLabel(slInfoTag, infoHandle);
00284 
00285 // Look first at the jetTags
00286 
00287   for (unsigned int iJetLabel = 0; iJetLabel != jetTagInputTags.size(); ++iJetLabel) {
00288     edm::Handle<reco::JetTagCollection> tagHandle;
00289     iEvent.getByLabel(jetTagInputTags[iJetLabel], tagHandle);
00290     const reco::JetTagCollection & tagColl = *(tagHandle.product());
00291     LogDebug("Info") << "Found " << tagColl.size() << " B candidates in collection " << jetTagInputTags[iJetLabel];
00292 
00293     int plotterSize =  binJetTagPlotters[iJetLabel].size();
00294     for (JetTagCollection::const_iterator tagI = tagColl.begin();
00295          tagI != tagColl.end(); ++tagI) {
00296       // Identify parton associated to jet.
00297 
00299       if (flavours[tagI->first] == 5 &&
00300           ((electronPlots && !leptons[tagI->first].electron) ||
00301            (muonPlots && !leptons[tagI->first].muon) ||
00302            (tauPlots && !leptons[tagI->first].tau)))
00303         continue;
00304 
00305       JetWithFlavour jetWithFlavour;
00306       if (!getJetWithFlavour(tagI->first, flavours, jetWithFlavour, iSetup))
00307         continue;
00308       if (!jetSelector(jetWithFlavour.first, std::abs(jetWithFlavour.second.getFlavour()), infoHandle))
00309         continue;
00310 
00311       for (int iPlotter = 0; iPlotter != plotterSize; ++iPlotter) {
00312               bool inBin = false;
00313               if (partonKinematics)
00314                 inBin = binJetTagPlotters[iJetLabel][iPlotter]->etaPtBin().inBin(jetWithFlavour.second.getLorentzVector().Eta(),
00315                                                                                  jetWithFlavour.second.getLorentzVector().Pt());
00316               else
00317                 inBin = binJetTagPlotters[iJetLabel][iPlotter]->etaPtBin().inBin(jetWithFlavour.first);
00318               // Fill histograms if in desired pt/rapidity bin.
00319               if (inBin)
00320                 binJetTagPlotters[iJetLabel][iPlotter]->analyzeTag(jetWithFlavour.first, tagI->second, std::abs(jetWithFlavour.second.getFlavour()));
00321       }
00322     }
00323   }
00324 
00325 // Now look at Tag Correlations
00326   for (unsigned int iJetLabel = 0; iJetLabel != tagCorrelationInputTags.size(); ++iJetLabel) {
00327     const std::pair<edm::InputTag, edm::InputTag>& inputTags = tagCorrelationInputTags[iJetLabel];
00328     edm::Handle<reco::JetTagCollection> tagHandle1;
00329     iEvent.getByLabel(inputTags.first, tagHandle1);
00330     const reco::JetTagCollection& tagColl1 = *(tagHandle1.product());
00331 
00332     edm::Handle<reco::JetTagCollection> tagHandle2;
00333     iEvent.getByLabel(inputTags.second, tagHandle2);
00334     const reco::JetTagCollection& tagColl2 = *(tagHandle2.product());
00335 
00336     int plotterSize = binTagCorrelationPlotters[iJetLabel].size();
00337     for (JetTagCollection::const_iterator tagI = tagColl1.begin(); tagI != tagColl1.end(); ++tagI) {
00338       
00339       if (flavours[tagI->first] == 5 &&
00340           ((electronPlots && !leptons[tagI->first].electron) ||
00341            (muonPlots && !leptons[tagI->first].muon) ||
00342            (tauPlots && !leptons[tagI->first].tau)))
00343         continue;
00344 
00345       JetWithFlavour jetWithFlavour;
00346       if (!getJetWithFlavour(tagI->first, flavours, jetWithFlavour, iSetup))
00347         continue;
00348       if (!jetSelector(jetWithFlavour.first, std::abs(jetWithFlavour.second.getFlavour()), infoHandle))
00349         continue;
00350 
00351       for(int iPlotter = 0; iPlotter != plotterSize; ++iPlotter) {
00352         bool inBin = false;
00353         if (partonKinematics)
00354           inBin = binTagCorrelationPlotters[iJetLabel][iPlotter]->etaPtBin().inBin(jetWithFlavour.second.getLorentzVector().Eta(),
00355                                                                                    jetWithFlavour.second.getLorentzVector().Pt());
00356 
00357         else
00358           inBin = binTagCorrelationPlotters[iJetLabel][iPlotter]->etaPtBin().inBin(jetWithFlavour.first);
00359 
00360         if(inBin)
00361         {
00362           double discr2 = tagColl2[tagI->first];
00363           binTagCorrelationPlotters[iJetLabel][iPlotter]->analyzeTags(tagI->second, discr2, std::abs(jetWithFlavour.second.getFlavour()));
00364         }
00365       }
00366     }
00367   }
00368 
00369 // Now look at the TagInfos
00370 
00371   for (unsigned int iJetLabel = 0; iJetLabel != tiDataFormatType.size(); ++iJetLabel) {
00372     int plotterSize = binTagInfoPlotters[iJetLabel].size();
00373     for (int iPlotter = 0; iPlotter != plotterSize; ++iPlotter)
00374       binTagInfoPlotters[iJetLabel][iPlotter]->setEventSetup(iSetup);
00375 
00376     vector<edm::InputTag> & inputTags = tagInfoInputTags[iJetLabel];
00377     if (inputTags.empty()) {
00378       // deferred retrieval of input tags
00379       BaseTagInfoPlotter *firstPlotter = binTagInfoPlotters[iJetLabel][0];
00380       int iModule = binTagInfoPlottersToModuleConfig[firstPlotter];
00381       vector<string> labels = firstPlotter->tagInfoRequirements();
00382       if (labels.empty())
00383         labels.push_back("label");
00384       for (vector<string>::const_iterator iLabels = labels.begin();
00385            iLabels != labels.end(); ++iLabels) {
00386         edm::InputTag inputTag =
00387                 moduleConfig[iModule].getParameter<InputTag>(*iLabels);
00388         inputTags.push_back(inputTag);
00389       }
00390     }
00391 
00392     unsigned int nInputTags = inputTags.size();
00393     vector< edm::Handle< View<BaseTagInfo> > > tagInfoHandles(nInputTags);
00394     edm::ProductID jetProductID;
00395     unsigned int nTagInfos = 0;
00396     for (unsigned int iInputTags = 0; iInputTags < inputTags.size(); ++iInputTags) {
00397       edm::Handle< View<BaseTagInfo> > & tagInfoHandle = tagInfoHandles[iInputTags];
00398       iEvent.getByLabel(inputTags[iInputTags], tagInfoHandle);
00399       unsigned int size = tagInfoHandle->size();
00400       LogDebug("Info") << "Found " << size << " B candidates in collection " << inputTags[iInputTags];
00401 
00402       edm::ProductID thisProductID = (size > 0) ? (*tagInfoHandle)[0].jet().id() : edm::ProductID();
00403       if (iInputTags == 0) {
00404         jetProductID = thisProductID;
00405         nTagInfos = size;
00406       } else if (jetProductID != thisProductID)
00407         throw cms::Exception("Configuration") << "TagInfos are referencing a different jet collection." << endl;
00408       else if (nTagInfos != size)
00409         throw cms::Exception("Configuration") << "TagInfo collections are having a different size." << endl;
00410     }
00411 
00412     for (unsigned int iTagInfos = 0; iTagInfos < nTagInfos; ++iTagInfos) {
00413       vector<const BaseTagInfo*> baseTagInfos(nInputTags);
00414       edm::RefToBase<Jet> jetRef;
00415       for (unsigned int iTagInfo = 0; iTagInfo < nInputTags; iTagInfo++) {
00416         const BaseTagInfo &baseTagInfo = (*tagInfoHandles[iTagInfo])[iTagInfos];
00417         if (iTagInfo == 0)
00418           jetRef = baseTagInfo.jet();
00419         else if (baseTagInfo.jet() != jetRef)
00420           throw cms::Exception("Configuration") << "TagInfos pointing to different jets." << endl;
00421         baseTagInfos[iTagInfo] = &baseTagInfo;
00422       }
00423 
00424       // Identify parton associated to jet.
00425 
00427       if (flavours[jetRef] == 5 &&
00428           ((electronPlots && !leptons[jetRef].electron) ||
00429            (muonPlots && !leptons[jetRef].muon) ||
00430            (tauPlots && !leptons[jetRef].tau)))
00431          continue;
00432 
00433       JetWithFlavour jetWithFlavour;
00434       if (!getJetWithFlavour(jetRef, flavours, jetWithFlavour, iSetup))
00435         continue;
00436       if (!jetSelector(jetWithFlavour.first, std::abs(jetWithFlavour.second.getFlavour()), infoHandle))
00437         continue;
00438 
00439       for (int iPlotter = 0; iPlotter != plotterSize; ++iPlotter) {
00440               bool inBin = false;
00441               if (partonKinematics)
00442                inBin = binTagInfoPlotters[iJetLabel][iPlotter]->etaPtBin().inBin(jetWithFlavour.second.getLorentzVector().Eta(),
00443                                                                                  jetWithFlavour.second.getLorentzVector().Pt());
00444               else
00445                inBin = binTagInfoPlotters[iJetLabel][iPlotter]->etaPtBin().inBin(*jetRef);
00446               // Fill histograms if in desired pt/rapidity bin.
00447               if (inBin)
00448                 binTagInfoPlotters[iJetLabel][iPlotter]->analyzeTag(baseTagInfos, std::abs(jetWithFlavour.second.getFlavour()));
00449       }
00450     }
00451   }
00452 }
00453 
00454 bool  BTagPerformanceAnalyzerMC::getJetWithFlavour(     edm::RefToBase<Jet> jetRef, FlavourMap flavours,
00455         JetWithFlavour & jetWithFlavour, const edm::EventSetup & es)
00456 {
00457 
00458   edm::ProductID recProdId = jetRef.id();
00459   edm::ProductID refProdId = (flavours.begin() == flavours.end())
00460     ? recProdId
00461     : flavours.begin()->first.id();
00462 
00463   if (!eventInitialized) {
00464     jetCorrector.setEventSetup(es);
00465     if (recProdId != refProdId) {
00466       edm::RefToBaseVector<Jet> refJets;
00467       for(FlavourMap::const_iterator iter = flavours.begin();
00468           iter != flavours.end(); ++iter)
00469         refJets.push_back(iter->first);
00470       const edm::RefToBaseProd<Jet> recJetsProd(jetRef);
00471       edm::RefToBaseVector<Jet> recJets;
00472       for(unsigned int i = 0; i < recJetsProd->size(); i++)
00473         recJets.push_back(edm::RefToBase<Jet>(recJetsProd, i));
00474       jetMatcher.matchCollections(refJets, recJets, es);
00475     }
00476     eventInitialized = true;
00477   }
00478 
00479   if (recProdId != refProdId) {
00480     jetRef = jetMatcher(jetRef);
00481     if (jetRef.isNull())
00482       return false;
00483   }
00484 
00485   jetWithFlavour.first = jetCorrector(*jetRef);
00486 
00487   jetWithFlavour.second = reco::JetFlavour(jetWithFlavour.first.p4(), math::XYZPoint (0,0,0), flavours[jetRef]);
00488 
00489   LogTrace("Info") << "Found jet with flavour "<<jetWithFlavour.second.getFlavour()<<endl;
00490   LogTrace("Info") << jetWithFlavour.first.p()<<" , "<< jetWithFlavour.first.pt()<<" - "
00491    << jetWithFlavour.second.getLorentzVector().P()<<" , "<< jetWithFlavour.second.getLorentzVector().Pt()<<endl;
00492 
00493   return true;
00494 }
00495 
00496 void BTagPerformanceAnalyzerMC::endJob()
00497 {
00498   if (!finalize) return;
00499   setTDRStyle();
00500   for (unsigned int iJetLabel = 0; iJetLabel != binJetTagPlotters.size(); ++iJetLabel) {
00501     int plotterSize =  binJetTagPlotters[iJetLabel].size();
00502     for (int iPlotter = 0; iPlotter != plotterSize; ++iPlotter) {
00503        binJetTagPlotters[iJetLabel][iPlotter]->finalize();
00504       //      binJetTagPlotters[iJetLabel][iPlotter]->write(allHisto);
00505       if (producePs)  (*binJetTagPlotters[iJetLabel][iPlotter]).psPlot(psBaseName);
00506       if (produceEps) (*binJetTagPlotters[iJetLabel][iPlotter]).epsPlot(epsBaseName);
00507     }
00508    
00509       if(makeDiffPlots_) { 
00510         for (vector<BTagDifferentialPlot *>::iterator iPlotter = differentialPlots[iJetLabel].begin();
00511              iPlotter != differentialPlots[iJetLabel].end(); ++ iPlotter) {
00512           (*iPlotter)->process();
00513           if (producePs)  (*iPlotter)->psPlot(psBaseName);
00514           if (produceEps) (*iPlotter)->epsPlot(epsBaseName);
00515           //      (**iPlotter).write(allHisto);
00516         }
00517       }
00518   }
00519   for (vector<vector<BaseTagInfoPlotter*> >::iterator iJetLabel = binTagInfoPlotters.begin();
00520        iJetLabel != binTagInfoPlotters.end(); ++iJetLabel) {
00521     for (vector<BaseTagInfoPlotter*>::iterator iPlotter = iJetLabel->begin(); iPlotter != iJetLabel->end(); ++iPlotter) {
00522       (*iPlotter)->finalize();
00523       //      binTagInfoPlotters[iJetLabel][iPlotter]->write(allHisto);
00524       if (producePs)  (*iPlotter)->psPlot(psBaseName);
00525       if (produceEps) (*iPlotter)->epsPlot(epsBaseName);
00526     }
00527   }
00528 }
00529 
00530 
00531 //define this as a plug-in
00532 DEFINE_FWK_MODULE(BTagPerformanceAnalyzerMC);