Go to the documentation of this file.00001 #include "TopQuarkAnalysis/TopJetCombination/plugins/TtSemiLepJetCombGeom.h"
00002
00003 #include "Math/VectorUtil.h"
00004
00005 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLepEvtPartons.h"
00006 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00008
00009 TtSemiLepJetCombGeom::TtSemiLepJetCombGeom(const edm::ParameterSet& cfg):
00010 jets_ (cfg.getParameter<edm::InputTag>("jets" )),
00011 leps_ (cfg.getParameter<edm::InputTag>("leps" )),
00012 maxNJets_ (cfg.getParameter<int> ("maxNJets" )),
00013 useDeltaR_ (cfg.getParameter<bool> ("useDeltaR" )),
00014 useBTagging_ (cfg.getParameter<bool> ("useBTagging" )),
00015 bTagAlgorithm_ (cfg.getParameter<std::string> ("bTagAlgorithm" )),
00016 minBDiscBJets_ (cfg.getParameter<double> ("minBDiscBJets" )),
00017 maxBDiscLightJets_(cfg.getParameter<double> ("maxBDiscLightJets"))
00018 {
00019 if(maxNJets_<4 && maxNJets_!=-1)
00020 throw cms::Exception("WrongConfig")
00021 << "Parameter maxNJets can not be set to " << maxNJets_ << ". \n"
00022 << "It has to be larger than 4 or can be set to -1 to take all jets.";
00023
00024 produces<std::vector<std::vector<int> > >();
00025 produces<int>("NumberOfConsideredJets");
00026 }
00027
00028 TtSemiLepJetCombGeom::~TtSemiLepJetCombGeom()
00029 {
00030 }
00031
00032 void
00033 TtSemiLepJetCombGeom::produce(edm::Event& evt, const edm::EventSetup& setup)
00034 {
00035 std::auto_ptr<std::vector<std::vector<int> > > pOut(new std::vector<std::vector<int> >);
00036 std::auto_ptr<int> pJetsConsidered(new int);
00037
00038 std::vector<int> match;
00039 for(unsigned int i = 0; i < 4; ++i)
00040 match.push_back( -1 );
00041
00042
00043 edm::Handle< std::vector<pat::Jet> > jets;
00044 evt.getByLabel(jets_, jets);
00045
00046
00047 edm::Handle< edm::View<reco::RecoCandidate> > leps;
00048 evt.getByLabel(leps_, leps);
00049
00050
00051 if(leps->empty() || jets->size() < 4){
00052 pOut->push_back( match );
00053 evt.put(pOut);
00054 *pJetsConsidered = jets->size();
00055 evt.put(pJetsConsidered, "NumberOfConsideredJets");
00056 return;
00057 }
00058
00059 unsigned maxNJets = maxNJets_;
00060 if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size();
00061 *pJetsConsidered = maxNJets;
00062 evt.put(pJetsConsidered, "NumberOfConsideredJets");
00063
00064 std::vector<bool> isBJet;
00065 std::vector<bool> isLJet;
00066 int cntBJets = 0;
00067 if(useBTagging_) {
00068 for(unsigned int idx=0; idx<maxNJets; ++idx) {
00069 isBJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ ) );
00070 isLJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_) );
00071 if((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ )cntBJets++;
00072 }
00073 }
00074
00075
00076
00077
00078
00079 double minDist=-1.;
00080 int lightQ =-1;
00081 int lightQBar=-1;
00082 for(unsigned int idx=0; idx<maxNJets; ++idx){
00083 if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue;
00084 for(unsigned int jdx=(idx+1); jdx<maxNJets; ++jdx){
00085 if(useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx]))) continue;
00086 double dist = distance((*jets)[idx].p4(), (*jets)[jdx].p4());
00087 if( minDist<0. || dist<minDist ){
00088 minDist=dist;
00089 lightQ =idx;
00090 lightQBar=jdx;
00091 }
00092 }
00093 }
00094
00095 reco::Particle::LorentzVector wHad;
00096 if( isValid(lightQ, jets) && isValid(lightQBar, jets) )
00097 wHad = (*jets)[lightQ].p4() + (*jets)[lightQBar].p4();
00098
00099
00100
00101
00102
00103 minDist=-1.;
00104 int hadB=-1;
00105 if( isValid(lightQ, jets) && isValid(lightQBar, jets) ) {
00106 for(unsigned int idx=0; idx<maxNJets; ++idx){
00107 if(useBTagging_ && !isBJet[idx]) continue;
00108
00109 if( (int)idx!=lightQ && (int)idx!=lightQBar ) {
00110 double dist = distance((*jets)[idx].p4(), wHad);
00111 if( minDist<0. || dist<minDist ){
00112 minDist=dist;
00113 hadB=idx;
00114 }
00115 }
00116 }
00117 }
00118
00119
00120
00121
00122
00123 minDist=-1.;
00124 int lepB=-1;
00125 for(unsigned int idx=0; idx<maxNJets; ++idx){
00126 if(useBTagging_ && !isBJet[idx]) continue;
00127
00128 if( (int)idx!=lightQ && (int)idx!=lightQBar && (int)idx!=hadB ){
00129 double dist = distance((*jets)[idx].p4(), (*leps)[0].p4());
00130 if( minDist<0. || dist<minDist ){
00131 minDist=dist;
00132 lepB=idx;
00133 }
00134 }
00135 }
00136
00137 match[TtSemiLepEvtPartons::LightQ ] = lightQ;
00138 match[TtSemiLepEvtPartons::LightQBar] = lightQBar;
00139 match[TtSemiLepEvtPartons::HadB ] = hadB;
00140 match[TtSemiLepEvtPartons::LepB ] = lepB;
00141
00142 pOut->push_back( match );
00143 evt.put(pOut);
00144 }
00145
00146 double
00147 TtSemiLepJetCombGeom::distance(const math::XYZTLorentzVector& v1, const math::XYZTLorentzVector& v2)
00148 {
00149
00150
00151 if(useDeltaR_) return ROOT::Math::VectorUtil::DeltaR(v1, v2);
00152 return fabs(v1.theta() - v2.theta());
00153 }