Go to the documentation of this file.00001 #include "TopQuarkAnalysis/TopJetCombination/plugins/TtSemiLepJetCombWMassDeltaTopMass.h"
00002
00003 #include "AnalysisDataFormats/TopObjects/interface/TtSemiLepEvtPartons.h"
00004 #include "DataFormats/PatCandidates/interface/Electron.h"
00005 #include "DataFormats/PatCandidates/interface/Muon.h"
00006 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00007 #include "TopQuarkAnalysis/TopTools/interface/MEzCalculator.h"
00008
00009 TtSemiLepJetCombWMassDeltaTopMass::TtSemiLepJetCombWMassDeltaTopMass(const edm::ParameterSet& cfg):
00010 jets_ (cfg.getParameter<edm::InputTag>("jets" )),
00011 leps_ (cfg.getParameter<edm::InputTag>("leps" )),
00012 mets_ (cfg.getParameter<edm::InputTag>("mets" )),
00013 maxNJets_ (cfg.getParameter<int> ("maxNJets" )),
00014 wMass_ (cfg.getParameter<double> ("wMass" )),
00015 useBTagging_ (cfg.getParameter<bool> ("useBTagging" )),
00016 bTagAlgorithm_ (cfg.getParameter<std::string> ("bTagAlgorithm" )),
00017 minBDiscBJets_ (cfg.getParameter<double> ("minBDiscBJets" )),
00018 maxBDiscLightJets_(cfg.getParameter<double> ("maxBDiscLightJets")),
00019 neutrinoSolutionType_(cfg.getParameter<int> ("neutrinoSolutionType"))
00020 {
00021 if(maxNJets_<4 && maxNJets_!=-1)
00022 throw cms::Exception("WrongConfig")
00023 << "Parameter maxNJets can not be set to " << maxNJets_ << ". \n"
00024 << "It has to be larger than 4 or can be set to -1 to take all jets.";
00025
00026 produces<std::vector<std::vector<int> > >();
00027 produces<int>("NumberOfConsideredJets");
00028 }
00029
00030 TtSemiLepJetCombWMassDeltaTopMass::~TtSemiLepJetCombWMassDeltaTopMass()
00031 {
00032 }
00033
00034 void
00035 TtSemiLepJetCombWMassDeltaTopMass::produce(edm::Event& evt, const edm::EventSetup& setup)
00036 {
00037 std::auto_ptr<std::vector<std::vector<int> > > pOut(new std::vector<std::vector<int> >);
00038 std::auto_ptr<int> pJetsConsidered(new int);
00039
00040 std::vector<int> match;
00041 for(unsigned int i = 0; i < 4; ++i)
00042 match.push_back( -1 );
00043
00044
00045 edm::Handle< std::vector<pat::Jet> > jets;
00046 evt.getByLabel(jets_, jets);
00047
00048
00049 edm::Handle< edm::View<reco::RecoCandidate> > leps;
00050 evt.getByLabel(leps_, leps);
00051
00052
00053 edm::Handle< std::vector<pat::MET> > mets;
00054 evt.getByLabel(mets_, mets);
00055
00056
00057 if(leps->empty() || jets->size() < 4 || mets->empty()){
00058 pOut->push_back( match );
00059 evt.put(pOut);
00060 *pJetsConsidered = jets->size();
00061 evt.put(pJetsConsidered, "NumberOfConsideredJets");
00062 return;
00063 }
00064
00065 unsigned maxNJets = maxNJets_;
00066 if(maxNJets_ == -1 || (int)jets->size() < maxNJets_) maxNJets = jets->size();
00067 *pJetsConsidered = maxNJets;
00068 evt.put(pJetsConsidered, "NumberOfConsideredJets");
00069
00070 std::vector<bool> isBJet;
00071 std::vector<bool> isLJet;
00072 int cntBJets = 0;
00073 if(useBTagging_) {
00074 for(unsigned int idx=0; idx<maxNJets; ++idx) {
00075 isBJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ ) );
00076 isLJet.push_back( ((*jets)[idx].bDiscriminator(bTagAlgorithm_) < maxBDiscLightJets_) );
00077 if((*jets)[idx].bDiscriminator(bTagAlgorithm_) > minBDiscBJets_ )cntBJets++;
00078 }
00079 }
00080
00081
00082
00083
00084
00085 double wDist =-1.;
00086 std::vector<int> closestToWMassIndices;
00087 closestToWMassIndices.push_back(-1);
00088 closestToWMassIndices.push_back(-1);
00089 for(unsigned idx=0; idx<maxNJets; ++idx){
00090 if(useBTagging_ && (!isLJet[idx] || (cntBJets<=2 && isBJet[idx]))) continue;
00091 for(unsigned jdx=(idx+1); jdx<maxNJets; ++jdx){
00092 if(useBTagging_ && (!isLJet[jdx] || (cntBJets<=2 && isBJet[jdx]) || (cntBJets==3 && isBJet[idx] && isBJet[jdx]))) continue;
00093 reco::Particle::LorentzVector sum =
00094 (*jets)[idx].p4()+
00095 (*jets)[jdx].p4();
00096 if( wDist<0. || wDist>fabs(sum.mass()-wMass_) ){
00097 wDist=fabs(sum.mass()-wMass_);
00098 closestToWMassIndices.clear();
00099 closestToWMassIndices.push_back(idx);
00100 closestToWMassIndices.push_back(jdx);
00101 }
00102 }
00103 }
00104
00105
00106
00107
00108 double neutrino_pz = 0.;
00109 if(neutrinoSolutionType_!=-1) {
00110 MEzCalculator mez;
00111 mez.SetMET( *mets->begin() );
00112 if( dynamic_cast<const reco::Muon*>(&(leps->front())) )
00113 mez.SetLepton( (*leps)[0], true );
00114 else if( dynamic_cast<const reco::GsfElectron*>(&(leps->front())) )
00115 mez.SetLepton( (*leps)[0], false );
00116 else
00117 throw cms::Exception("UnimplementedFeature") << "Type of lepton given together with MET for solving neutrino kinematics is neither muon nor electron.\n";
00118 neutrino_pz = mez.Calculate(neutrinoSolutionType_);
00119 }
00120 const math::XYZTLorentzVector neutrino( mets->begin()->px(), mets->begin()->py(), neutrino_pz, sqrt(mets->begin()->px()*mets->begin()->px()
00121 + mets->begin()->py()*mets->begin()->py()
00122 + neutrino_pz*neutrino_pz) );
00123 const reco::Particle::LorentzVector lepW = neutrino + leps->front().p4();
00124
00125
00126
00127
00128
00129
00130 double deltaTop=-1.;
00131 int hadB=-1;
00132 int lepB=-1;
00133 if( isValid(closestToWMassIndices[0], jets) && isValid(closestToWMassIndices[1], jets)) {
00134 const reco::Particle::LorentzVector hadW =
00135 (*jets)[closestToWMassIndices[0]].p4()+
00136 (*jets)[closestToWMassIndices[1]].p4();
00137
00138 for(unsigned idx=0; idx<maxNJets; ++idx){
00139 if(useBTagging_ && !isBJet[idx]) continue;
00140
00141 if( (int)idx!=closestToWMassIndices[0] && (int)idx!=closestToWMassIndices[1] ){
00142 reco::Particle::LorentzVector hadTop = hadW + (*jets)[idx].p4();
00143
00144 for(unsigned jdx=0; jdx<maxNJets; ++jdx){
00145 if(useBTagging_ && !isBJet[jdx]) continue;
00146
00147 if( (int)jdx!=closestToWMassIndices[0] && (int)jdx!=closestToWMassIndices[1] && jdx!=idx ){
00148 reco::Particle::LorentzVector lepTop = lepW + (*jets)[jdx].p4();
00149 if( deltaTop<0. || deltaTop>fabs(hadTop.mass()-lepTop.mass()) ){
00150 deltaTop=fabs(hadTop.mass()-lepTop.mass());
00151 hadB=idx;
00152 lepB=jdx;
00153 }
00154 }
00155 }
00156 }
00157 }
00158 }
00159
00160 match[TtSemiLepEvtPartons::LightQ ] = closestToWMassIndices[0];
00161 match[TtSemiLepEvtPartons::LightQBar] = closestToWMassIndices[1];
00162 match[TtSemiLepEvtPartons::HadB ] = hadB;
00163 match[TtSemiLepEvtPartons::LepB ] = lepB;
00164
00165 pOut->push_back( match );
00166 evt.put(pOut);
00167 }