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 }
00026
00027 TtSemiLepJetCombGeom::~TtSemiLepJetCombGeom()
00028 {
00029 }
00030
00031 void
00032 TtSemiLepJetCombGeom::produce(edm::Event& evt, const edm::EventSetup& setup)
00033 {
00034 std::auto_ptr< std::vector<std::vector<int> > >pOut (new std::vector<std::vector<int> >);
00035
00036 std::vector<int> match;
00037 for(unsigned int i = 0; i < 4; ++i)
00038 match.push_back( -1 );
00039
00040
00041 edm::Handle< std::vector<pat::Jet> > jets;
00042 evt.getByLabel(jets_, jets);
00043
00044
00045 edm::Handle< edm::View<reco::RecoCandidate> > leps;
00046 evt.getByLabel(leps_, leps);
00047
00048
00049 if(leps->empty() || jets->size() < 4){
00050 pOut->push_back( match );
00051 evt.put(pOut);
00052 return;
00053 }
00054
00055 unsigned maxNJets = maxNJets_;
00056 if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size();
00057
00058 std::vector<bool> isBJet;
00059 std::vector<bool> isLJet;
00060 int cntBJets = 0;
00061 if(useBTagging_) {
00062 for(unsigned int idx=0; idx<maxNJets; ++idx) {
00063 isBJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ ) );
00064 isLJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_) );
00065 if((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ )cntBJets++;
00066 }
00067 }
00068
00069
00070
00071
00072
00073 double minDist=-1.;
00074 int lightQ =-1;
00075 int lightQBar=-1;
00076 for(unsigned int idx=0; idx<maxNJets; ++idx){
00077 if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue;
00078 for(unsigned int jdx=(idx+1); jdx<maxNJets; ++jdx){
00079 if(useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx]))) continue;
00080 double dist = distance((*jets)[idx].p4(), (*jets)[jdx].p4());
00081 if( minDist<0. || dist<minDist ){
00082 minDist=dist;
00083 lightQ =idx;
00084 lightQBar=jdx;
00085 }
00086 }
00087 }
00088
00089 reco::Particle::LorentzVector wHad;
00090 if( isValid(lightQ, jets) && isValid(lightQBar, jets) )
00091 wHad = (*jets)[lightQ].p4() + (*jets)[lightQBar].p4();
00092
00093
00094
00095
00096
00097 minDist=-1.;
00098 int hadB=-1;
00099 if( isValid(lightQ, jets) && isValid(lightQBar, jets) ) {
00100 for(unsigned int idx=0; idx<maxNJets; ++idx){
00101 if(useBTagging_ && !isBJet[idx]) continue;
00102
00103 if( (int)idx!=lightQ && (int)idx!=lightQBar ) {
00104 double dist = distance((*jets)[idx].p4(), wHad);
00105 if( minDist<0. || dist<minDist ){
00106 minDist=dist;
00107 hadB=idx;
00108 }
00109 }
00110 }
00111 }
00112
00113
00114
00115
00116
00117 minDist=-1.;
00118 int lepB=-1;
00119 for(unsigned int idx=0; idx<maxNJets; ++idx){
00120 if(useBTagging_ && !isBJet[idx]) continue;
00121
00122 if( (int)idx!=lightQ && (int)idx!=lightQBar && (int)idx!=hadB ){
00123 double dist = distance((*jets)[idx].p4(), (*leps)[0].p4());
00124 if( minDist<0. || dist<minDist ){
00125 minDist=dist;
00126 lepB=idx;
00127 }
00128 }
00129 }
00130
00131 match[TtSemiLepEvtPartons::LightQ ] = lightQ;
00132 match[TtSemiLepEvtPartons::LightQBar] = lightQBar;
00133 match[TtSemiLepEvtPartons::HadB ] = hadB;
00134 match[TtSemiLepEvtPartons::LepB ] = lepB;
00135
00136 pOut->push_back( match );
00137 evt.put(pOut);
00138 }
00139
00140 double
00141 TtSemiLepJetCombGeom::distance(const math::XYZTLorentzVector& v1, const math::XYZTLorentzVector& v2)
00142 {
00143
00144
00145 if(useDeltaR_) return ROOT::Math::VectorUtil::DeltaR(v1, v2);
00146 return fabs(v1.theta() - v2.theta());
00147 }