CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
StEvtSolutionMaker.cc
Go to the documentation of this file.
1 //
2 //
3 
4 #include <memory>
5 
7 
9 {
10  // configurables
11  electronSrcToken_ = mayConsume<std::vector<pat::Electron> >(iConfig.getParameter<edm::InputTag>("electronSource"));
12  muonSrcToken_ = mayConsume<std::vector<pat::Muon> >(iConfig.getParameter<edm::InputTag>("muonSource"));
13  metSrcToken_ = consumes<std::vector<pat::MET> >(iConfig.getParameter<edm::InputTag>("metSource"));
14  jetSrcToken_ = consumes<std::vector<pat::Jet> >(iConfig.getParameter<edm::InputTag>("jetSource"));
15  genEvtSrcToken_ = mayConsume<StGenEvent>(edm::InputTag("genEvt"));
16  leptonFlavour_ = iConfig.getParameter< std::string >("leptonFlavour");
17  jetCorrScheme_ = iConfig.getParameter<int> ("jetCorrectionScheme");
18  //jetInput_ = iConfig.getParameter< std::string > ("jetInput");
19  doKinFit_ = iConfig.getParameter< bool > ("doKinFit");
20  addLRJetComb_ = iConfig.getParameter< bool > ("addLRJetComb");
21  maxNrIter_ = iConfig.getParameter< int > ("maxNrIter");
22  maxDeltaS_ = iConfig.getParameter< double > ("maxDeltaS");
23  maxF_ = iConfig.getParameter< double > ("maxF");
24  jetParam_ = iConfig.getParameter<int> ("jetParametrisation");
25  lepParam_ = iConfig.getParameter<int> ("lepParametrisation");
26  metParam_ = iConfig.getParameter<int> ("metParametrisation");
27  constraints_ = iConfig.getParameter< std::vector<int> > ("constraints");
28  matchToGenEvt_ = iConfig.getParameter< bool > ("matchToGenEvt");
29 
30  // define kinfitter
31  if(doKinFit_){
32  myKinFitter = new StKinFitter(jetParam_, lepParam_, metParam_, maxNrIter_, maxDeltaS_, maxF_, constraints_);
33  }
34  // define what will be produced
35  produces<std::vector<StEvtSolution> >();
36 }
37 
39 {
40  if (doKinFit_) delete myKinFitter;
41 }
42 
44 {
45  //
46  // TopObject Selection
47  //
48 
49  // select lepton (the TtLepton vectors are, for the moment, sorted on pT)
50  bool leptonFound = false;
52  if(leptonFlavour_ == "muon"){
53  iEvent.getByToken(muonSrcToken_, muons);
54  if( muons->size() > 0 ) leptonFound = true;
55  }
57  if(leptonFlavour_ == "electron"){
58  iEvent.getByToken(electronSrcToken_, electrons);
59  if( electrons->size() > 0 ) leptonFound = true;
60  }
61 
62  // select MET (TopMET vector is sorted on ET)
63  bool metFound = false;
65  iEvent.getByToken(metSrcToken_, mets);
66  if( mets->size() > 0 ) metFound = true;
67 
68  // select Jets
69  bool jetsFound = false;
71  iEvent.getByToken(jetSrcToken_, jets);
72  unsigned int maxJets=2;//this has to become a custom-defined parameter (we may want 2 or 3 jets)
73  if (jets->size() >= 2) jetsFound = true;
74 
75  std::vector<StEvtSolution> *evtsols = new std::vector<StEvtSolution>();
76  if(leptonFound && metFound && jetsFound){
77  std::cout<<"constructing solutions"<<std::endl;
78  for (unsigned int b=0; b<maxJets; b++) {
79  for (unsigned int l=0; l<maxJets; l++) {
80  if(b!=l){ // to avoid double counting
81  StEvtSolution asol;
83  if(leptonFlavour_ == "muon") asol.setMuon(muons, 0);
84  if(leptonFlavour_ == "electron") asol.setElectron(electrons, 0);
85  asol.setNeutrino(mets, 0);
86  asol.setBottom(jets, b);
87  asol.setLight(jets, l);
88 
89  if(doKinFit_) asol = myKinFitter->addKinFitInfo(&asol);
90 
91  /* to be adapted to ST (Andrea)
92  if(addLRJetComb_){
93  asol.setPtrueCombExist(jetCombProbs[m].getPTrueCombExist(&afitsol));
94  asol.setPtrueBJetSel(jetCombProbs[m].getPTrueBJetSel(&afitsol));
95  asol.setPtrueBhadrSel(jetCombProbs[m].getPTrueBhadrSel(&afitsol));
96  asol.setPtrueJetComb(afitsol.getPtrueCombExist()*afitsol.getPtrueBJetSel()*afitsol.getPtrueBhadrSel());
97  }
98  */
99  evtsols->push_back(asol);
100  }
101  }
102  }
103 
104  // if asked for, match the event solutions to the gen Event
105  if(matchToGenEvt_){
106  /*
107  edm::Handle<StGenEvent> genEvt;
108  iEvent.getByToken (genEvtSrcToken_,genEvt);
109  double bestSolDR = 9999.;
110  int bestSol = 0;
111  for(size_t s=0; s<evtsols->size(); s++) {
112  (*evtsols)[s].setGenEvt(genEvt->particles());
113  vector<double> bm = BestMatch((*evtsols)[s], false); //false to use DR, true SpaceAngles
114  (*evtsols)[s].setSumDeltaRjp(bm[0]); // dRBB + dRLL
115  (*evtsols)[s].setChangeBL((int) bm[1]); // has swapped or not
116  (*evtsols)[s].setDeltaRB(bm[2]);
117  (*evtsols)[s].setDeltaRL(bm[3]);
118  if(bm[0]<bestSolDR) { bestSolDR = bm[0]; bestSol = s; }
119  }
120  (*evtsols)[bestSol].setBestSol(true);
121  */
122  }
123 
124  //store the vector of solutions to the event
125  std::auto_ptr<std::vector<StEvtSolution> > pOut(evtsols);
126  iEvent.put(pOut);
127  }
128  else
129  {
130 
131  std::cout<<"@@@ No calibrated solutions built, because: " << std::endl;;
132  if(jets->size()<maxJets) std::cout<<"@ nr jets = " << jets->size() << " < " << maxJets <<std::endl;
133  if(leptonFlavour_ == "muon" && !leptonFound) std::cout<<"@ no good muon candidate"<<std::endl;
134  if(leptonFlavour_ == "electron" && !leptonFound) std::cout<<"@ no good electron candidate"<<std::endl;
135  if(mets->size() == 0) std::cout<<"@ no MET reconstruction"<<std::endl;
136 
137  StEvtSolution asol;
138  evtsols->push_back(asol);
139  std::auto_ptr<std::vector<StEvtSolution> > pOut(evtsols);
140  iEvent.put(pOut);
141  }
142 }
T getParameter(std::string const &) const
void setMuon(const edm::Handle< std::vector< pat::Muon > > &muon, int i)
StKinFitter * myKinFitter
void setLight(const edm::Handle< std::vector< pat::Jet > > &jet, int i)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
void setNeutrino(const edm::Handle< std::vector< pat::MET > > &met, int i)
std::string leptonFlavour_
edm::EDGetTokenT< std::vector< pat::Jet > > jetSrcToken_
virtual void produce(edm::Event &, const edm::EventSetup &)
std::vector< int > constraints_
edm::EDGetTokenT< std::vector< pat::Muon > > muonSrcToken_
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
vector< PseudoJet > jets
edm::EDGetTokenT< StGenEvent > genEvtSrcToken_
StEvtSolutionMaker(const edm::ParameterSet &)
edm::EDGetTokenT< std::vector< pat::Electron > > electronSrcToken_
edm::EDGetTokenT< std::vector< pat::MET > > metSrcToken_
void setJetCorrectionScheme(int scheme)
double b
Definition: hdecay.h:120
void setBottom(const edm::Handle< std::vector< pat::Jet > > &jet, int i)
tuple muons
Definition: patZpeak.py:38
StEvtSolution addKinFitInfo(StEvtSolution *asol)
Definition: StKinFitter.cc:64
tuple cout
Definition: gather_cfg.py:121
void setElectron(const edm::Handle< std::vector< pat::Electron > > &elec, int i)