CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PatBTagCommonHistos.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PatBTag
4 // Class: PatBTagCommonHistos
5 //
14 //
15 // Original Author: J. E. Ramirez
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
28 
31 
32 //
33 // constructors and destructor
34 //
36  histocontainer_()
37 ,BTagger(iConfig.getParameter< edm::ParameterSet >("BJetOperatingPoints"))
38  //now do what ever initialization is needed
39 {
40  edm::ParameterSet PatBjet_(iConfig.getParameter< edm::ParameterSet >("BjetTag"));
41  BTagdisccut_ = PatBjet_.getUntrackedParameter<double>("mindiscriminatorcut",5.0);
42  BTagdiscriminator_= PatBjet_.getParameter<std::string>("discriminator");
43  BTagpurity_ = PatBjet_.getParameter<std::string>("purity");
44 }
45 
46 
48 {
49 
50  // do anything here that needs to be done at desctruction time
51  // (e.g. close files, deallocate resources etc.)
52 
53 }
54 
55 
56 //
57 // member functions
58 //
59 
60 // ------------ method called to for each event ------------
61 void
63 {
64 
65 float isb =jet_iter->bDiscriminator(BTagdiscriminator_);
66 //
67 // no pt cut histograms
68 //
69 if(
70  fabs(jet_iter->eta()) < 2.4
71  && fabs(jet_iter->eta()) > 0
72 ){
73  //
74  //Fill histos for Loose(defined in configuration file) cut
75  //tagged using auxiliar function (Loose)
76  //
78  histocontainer_["jet_pt_uncorr_"+flavor+"_tagged"]->Fill(jet_iter->correctedJet("raw").pt());
79  histocontainer_["jet_pt_"+flavor+"_tagged"]->Fill(jet_iter->pt());
80  histocontainer_["jet_eta_"+flavor+"_tagged"]->Fill(jet_iter->eta());
81  }
82  //
83  // Fill histos
84  // tagged minimum discriminant cut criteria
85  //
86  if( isb > float(BTagdisccut_) ){
87  histocontainer_["jet_pt_uncorr_"+flavor+"_taggedmincut"]->Fill(jet_iter->correctedJet("raw").pt());
88  histocontainer_["jet_pt_"+flavor+"_taggedmincut"]->Fill(jet_iter->pt());
89  histocontainer_["jet_eta_"+flavor+"_taggedmincut"]->Fill(jet_iter->eta());
90  }
91  //
92  //fill taggable jets histograms (tracks in jet > 0,1,2)
93  std::map<int,std::string> tagbl;
94  tagbl[0]="0";tagbl[1]="1";tagbl[2]="2";
95  for (size_t i=0; i < tagbl.size(); i++)
96  if( jet_iter->associatedTracks().size() > i ){
97  histocontainer_["jet_pt_uncorr_"+flavor+"_taggable"+tagbl[i]]->Fill(jet_iter->correctedJet("raw").pt());
98  histocontainer_["jet_pt_"+flavor+"_taggable"+tagbl[i]]->Fill(jet_iter->pt());
99  histocontainer_["jet_eta_"+flavor+"_taggable"+tagbl[i]]->Fill(jet_iter->eta());
100  if (jet_iter->pt() < 30 ){
101  histocontainer_["jet_eta_"+flavor+"_taggable"+tagbl[i]+"_030"]->Fill(jet_iter->eta());
102  }
103  if (jet_iter->pt() > 30 && jet_iter->pt() < 50){
104  histocontainer_["jet_eta_"+flavor+"_taggable"+tagbl[i]+"_3050"]->Fill(jet_iter->eta());
105  }
106  if (jet_iter->pt() > 50 ){
107  histocontainer_["jet_eta_"+flavor+"_taggable"+tagbl[i]+"_50"]->Fill(jet_iter->eta());
108  }
109  }
110  //
111  // Fill histos needed for normalization
112  // uncorrected pt distributions
113  // corrected pt distributions
114  // eta distributions
115  // scatter plots for pts
116  histocontainer_["jet_pt_uncorr_"+flavor]->Fill(jet_iter->correctedJet("raw").pt());
117  histocontainer_["jet_pt_"+flavor]->Fill(jet_iter->pt());
118  histocontainer_["jet_eta_"+flavor]->Fill(jet_iter->eta());
119  histocontainer_["tracks_in_jet_"+flavor]->Fill(jet_iter->associatedTracks().size());
120  h2_["jet_scatter_pt_"+flavor]->Fill(jet_iter->pt(),jet_iter->correctedJet("raw").pt());
121  for(size_t itrack=0;itrack < jet_iter->associatedTracks().size() ;++itrack){
122  histocontainer_["pt_tracks_in_jet_"+flavor]->Fill(jet_iter->associatedTracks()[itrack]->pt());
123  }
124  if (jet_iter->pt() < 30 ){
125  histocontainer_["jet_eta_"+flavor+"030"]->Fill(jet_iter->eta());
126  histocontainer_["tracks_in_jet_"+flavor+"030"]->Fill(jet_iter->associatedTracks().size());
127  for(size_t itrack=0;itrack < jet_iter->associatedTracks().size() ;++itrack){
128  histocontainer_["pt_tracks_in_jet_"+flavor+"030"]->Fill(jet_iter->associatedTracks()[itrack]->pt());
129  }
130  if (fabs(jet_iter->eta()) < 1.5 ){
131  for(size_t itrack=0;itrack < jet_iter->associatedTracks().size() ;++itrack){
132  histocontainer_["pt_tracks_in_jet_"+flavor+"_center030"]->Fill(jet_iter->associatedTracks()[itrack]->pt());
133  }
134  }
135  if (fabs(jet_iter->eta()) > 1.5 && fabs(jet_iter->eta()) < 2.4){
136  for(size_t itrack=0;itrack < jet_iter->associatedTracks().size() ;++itrack){
137  histocontainer_["pt_tracks_in_jet_"+flavor+"_sides030"]->Fill(jet_iter->associatedTracks()[itrack]->pt());
138  }
139  }
140  }
141  if (jet_iter->pt() > 30 && jet_iter->pt() < 50){
142  histocontainer_["jet_eta_"+flavor+"3050"]->Fill(jet_iter->eta());
143  histocontainer_["tracks_in_jet_"+flavor+"3050"]->Fill(jet_iter->associatedTracks().size());
144  for(size_t itrack=0;itrack < jet_iter->associatedTracks().size() ;++itrack){
145  histocontainer_["pt_tracks_in_jet_"+flavor+"3050"]->Fill(jet_iter->associatedTracks()[itrack]->pt());
146  }
147  if (fabs(jet_iter->eta()) < 1.5 ){
148  for(size_t itrack=0;itrack < jet_iter->associatedTracks().size() ;++itrack){
149  histocontainer_["pt_tracks_in_jet_"+flavor+"_center3050"]->Fill(jet_iter->associatedTracks()[itrack]->pt());
150  }
151  }
152  if (fabs(jet_iter->eta()) > 1.5 && fabs(jet_iter->eta()) < 2.4){
153  for(size_t itrack=0;itrack < jet_iter->associatedTracks().size() ;++itrack){
154  histocontainer_["pt_tracks_in_jet_"+flavor+"_sides3050"]->Fill(jet_iter->associatedTracks()[itrack]->pt());
155  }
156  }
157  }
158  if (jet_iter->pt() > 50 ){
159  histocontainer_["jet_eta_"+flavor+"50"]->Fill(jet_iter->eta());
160  histocontainer_["tracks_in_jet_"+flavor+"50"]->Fill(jet_iter->associatedTracks().size());
161  for(size_t itrack=0;itrack < jet_iter->associatedTracks().size() ;++itrack){
162  histocontainer_["pt_tracks_in_jet_"+flavor+"50"]->Fill(jet_iter->associatedTracks()[itrack]->pt());
163  }
164  if (fabs(jet_iter->eta()) < 1.5 ){
165  for(size_t itrack=0;itrack < jet_iter->associatedTracks().size() ;++itrack){
166  histocontainer_["pt_tracks_in_jet_"+flavor+"_center50"]->Fill(jet_iter->associatedTracks()[itrack]->pt());
167  }
168  }
169  if (fabs(jet_iter->eta()) > 1.5 && fabs(jet_iter->eta()) < 2.4){
170  for(size_t itrack=0;itrack < jet_iter->associatedTracks().size() ;++itrack){
171  histocontainer_["pt_tracks_in_jet_"+flavor+"_sides50"]->Fill(jet_iter->associatedTracks()[itrack]->pt());
172  }
173  }
174  }
175  if (fabs(jet_iter->eta()) < 1.5 ){
176  for(size_t itrack=0;itrack < jet_iter->associatedTracks().size() ;++itrack){
177  histocontainer_["pt_tracks_in_jet_"+flavor+"_center"]->Fill(jet_iter->associatedTracks()[itrack]->pt());
178  }
179  }
180  if (fabs(jet_iter->eta()) > 1.5 && fabs(jet_iter->eta()) < 2.4){
181  for(size_t itrack=0;itrack < jet_iter->associatedTracks().size() ;++itrack){
182  histocontainer_["pt_tracks_in_jet_"+flavor+"_sides"]->Fill(jet_iter->associatedTracks()[itrack]->pt());
183  }
184  }
185 
186  }//endif
187 }//end function
188 // ------------ method called once each job just before starting event loop ------------
189 // ------------ This function is needed to set a group of histogram -------------------
190 void
192 {
193 
194  const int ntptarray = 23;
195  const int njptarray = 14;
196  const int netaarray = 19;
197  Double_t jetetabins[netaarray] = {-2.5,-2.0,-1.75,-1.5,-1.25,-1.0,-0.75,-0.5,-0.25,0.0,0.25,0.5,0.75,1.0,1.25,1.5,1.75,2.0,2.5};
198  Double_t jetptxbins[njptarray] = {0.,10.,20.,30., 40., 50., 60., 70., 80, 90., 100., 120., 140., 230.};
199  Double_t jetpttbins[ntptarray] = {0.,1.,2.,4.,6.,8.,10.,15.,20.,25.,30., 35.,40., 45., 50., 60., 70., 80, 90., 100., 120., 140., 230.};
201  std::string histoid;
202  std::string histotitle;
203 
204  //Define histograms
205  histoid = "jet_pt_"+flavor; histotitle = flavor+" jet p_{T} [GeV/c]";
206  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),njptarray-1,jetptxbins);
207  histoid = "jet_pt_uncorr_"+flavor; histotitle = flavor+" jet uncorrected p_{T} [GeV/c]";
208  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),njptarray-1,jetptxbins);
209  histoid = "jet_eta_"+flavor;
210  histocontainer_["jet_eta_"+flavor]=fs->make<TH1D>(histoid.c_str(),"jet #eta",netaarray-1,jetetabins);
211  histoid = "jet_eta_"+flavor+"030";
212  histocontainer_["jet_eta_"+flavor+"030"]=fs->make<TH1D>(histoid.c_str(),"jet #eta",netaarray-1,jetetabins);
213  histoid = "jet_eta_"+flavor+"3050";
214  histocontainer_["jet_eta_"+flavor+"3050"]=fs->make<TH1D>(histoid.c_str(),"jet #eta",netaarray-1,jetetabins);
215  histoid = "jet_eta_"+flavor+"50";
216  histocontainer_["jet_eta_"+flavor+"50"]=fs->make<TH1D>(histoid.c_str(),"jet #eta",netaarray-1,jetetabins);
217  histoid = "jet_pt_"+flavor+"_tagged";
218  histocontainer_["jet_pt_"+flavor+"_tagged"]=fs->make<TH1D>(histoid.c_str(),"jet_tagged p_{T} [GeV/c]",njptarray-1,jetptxbins);
219  histoid = "jet_pt_uncorr_"+flavor+"_tagged";
220  histocontainer_["jet_pt_uncorr_"+flavor+"_tagged"]=fs->make<TH1D>(histoid.c_str(),"jet_tagged p_{T} [GeV/c]",njptarray-1,jetptxbins);
221  histoid = "jet_eta_"+flavor+"_tagged";
222  histocontainer_["jet_eta_"+flavor+"_tagged"]=fs->make<TH1D>(histoid.c_str(),"jet_tagged #eta",netaarray-1,jetetabins);
223 
224  //tagged min cut
225  histoid = "jet_pt_"+flavor+"_taggedmincut";
226  histocontainer_["jet_pt_"+flavor+"_taggedmincut"]=fs->make<TH1D>(histoid.c_str(),"jet_tagged p_{T} [GeV/c]",njptarray-1,jetptxbins);
227  histoid = "jet_pt_uncorr_"+flavor+"_taggedmincut";
228  histocontainer_["jet_pt_uncorr_"+flavor+"_taggedmincut"]=fs->make<TH1D>(histoid.c_str(),"jet_tagged p_{T} [GeV/c]",njptarray-1,jetptxbins);
229  histoid = "jet_eta_"+flavor+"_taggedmincut";
230  histocontainer_["jet_eta_"+flavor+"_taggedmincut"]=fs->make<TH1D>(histoid.c_str(),"jet_tagged #eta",netaarray-1,jetetabins);
231 
232  //#tracks per jet
233  histoid = "tracks_in_jet_"+flavor; histotitle = "traks per jet "+flavor;
234  histocontainer_["tracks_in_jet_"+flavor]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),31,-0.5,30.5);
235  histoid = "tracks_in_jet_"+flavor+"030"; histotitle = "traks per jet "+flavor+ "pt_{T} < 30[GeV/c]";
236  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),31,-0.5,30.5);
237  histoid = "tracks_in_jet_"+flavor+"3050"; histotitle = "traks per jet "+flavor+ "30 < pt_{T} < 50[GeV/c]";
238  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),31,-0.5,30.5);
239  histoid = "tracks_in_jet_"+flavor+"50"; histotitle = "traks per jet "+flavor+ "pt_{T} > 50[GeV/c]";
240  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),31,-0.5,30.5);
241 
242  // pt of tracks in bins of jet pt 0-30,30-50,50
243  histoid= "pt_tracks_in_jet_"+flavor; histotitle = "track p_{T} [GeV/c] "+ flavor+" jets";
244  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),ntptarray-1,jetpttbins);
245  histoid= "pt_tracks_in_jet_"+flavor+"030"; histotitle = "track p_{T} [GeV/c] "+ flavor+" jets";
246  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),ntptarray-1,jetpttbins);
247  histoid= "pt_tracks_in_jet_"+flavor+"3050"; histotitle = "track p_{T} [GeV/c] "+ flavor+" jets";
248  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),ntptarray-1,jetpttbins);
249  histoid= "pt_tracks_in_jet_"+flavor+"50"; histotitle = "track p_{T} [GeV/c] "+ flavor+" jets";
250  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),ntptarray-1,jetpttbins);
251 
252  // pt of tracks in bins of jet eta abs(eta)<1.5, 1.5<abs(eta)<2.4
253  // combined with bins of jet pt
254  histoid= "pt_tracks_in_jet_"+flavor+"_center"; histotitle = "track p_{T} [GeV/c] "+ flavor+" jets";
255  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),ntptarray-1,jetpttbins);
256  histoid= "pt_tracks_in_jet_"+flavor+"_center030"; histotitle = "track p_{T} [GeV/c] "+ flavor+" jets";
257  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),ntptarray-1,jetpttbins);
258  histoid= "pt_tracks_in_jet_"+flavor+"_center3050"; histotitle = "track p_{T} [GeV/c] "+ flavor+" jets";
259  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),ntptarray-1,jetpttbins);
260  histoid= "pt_tracks_in_jet_"+flavor+"_center50"; histotitle = "track p_{T} [GeV/c] "+ flavor+" jets";
261  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),ntptarray-1,jetpttbins);
262  histoid= "pt_tracks_in_jet_"+flavor+"_sides"; histotitle = "track p_{T} [GeV/c] "+ flavor+" jets";
263  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),ntptarray-1,jetpttbins);
264  histoid= "pt_tracks_in_jet_"+flavor+"_sides030"; histotitle = "track p_{T} [GeV/c] "+ flavor+" jets";
265  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),ntptarray-1,jetpttbins);
266  histoid= "pt_tracks_in_jet_"+flavor+"_sides3050"; histotitle = "track p_{T} [GeV/c] "+ flavor+" jets";
267  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),ntptarray-1,jetpttbins);
268  histoid= "pt_tracks_in_jet_"+flavor+"_sides50"; histotitle = "track p_{T} [GeV/c] "+ flavor+" jets";
269  histocontainer_[histoid]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),ntptarray-1,jetpttbins);
270 
271  //taggable 0,1,2
272  std::map<int,std::string> tagbl;
273  tagbl[0]="0";tagbl[1]="1";tagbl[2]="2";
274  for (size_t i=0; i < tagbl.size(); i++){
275  histoid = "jet_pt_"+flavor+"_taggable"+tagbl[i]; histotitle = flavor+" jet_taggable p_{T} [GeV/c]";
276  histocontainer_["jet_pt_"+flavor+"_taggable"+tagbl[i]]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),njptarray-1,jetptxbins);
277  histoid = "jet_pt_uncorr_"+flavor+"_taggable"+tagbl[i]; histotitle = flavor+" jet_taggable p_{T} uncorr [GeV/c]";
278  histocontainer_["jet_pt_uncorr_"+flavor+"_taggable"+tagbl[i]]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),njptarray-1,jetptxbins);
279  histoid = "jet_eta_"+flavor+"_taggable"+tagbl[i]; histotitle = flavor+" jet_taggable #eta";
280  histocontainer_["jet_eta_"+flavor+"_taggable"+tagbl[i]]=fs->make<TH1D>(histoid.c_str(),histotitle.c_str(),netaarray-1,jetetabins);
281  histoid = "jet_eta_"+flavor+"_taggable"+tagbl[i]+"_030";
282  histocontainer_["jet_eta_"+flavor+"_taggable"+tagbl[i]+"_030"]=fs->make<TH1D>(histoid.c_str(),"jet #eta",netaarray-1,jetetabins);
283  histoid = "jet_eta_"+flavor+"_taggable"+tagbl[i]+"_3050";
284  histocontainer_["jet_eta_"+flavor+"_taggable"+tagbl[i]+"_3050"]=fs->make<TH1D>(histoid.c_str(),"jet #eta",netaarray-1,jetetabins);
285  histoid = "jet_eta_"+flavor+"_taggable"+tagbl[i]+"_50";
286  histocontainer_["jet_eta_"+flavor+"_taggable"+tagbl[i]+"_50"]=fs->make<TH1D>(histoid.c_str(),"jet #eta",netaarray-1,jetetabins);
287  }
288 
289  histoid = "jet_scatter_pt_"+flavor;
290  h2_["jet_scatter_pt_"+flavor]=fs->make<TH2D>(histoid.c_str(),"jet p_{T} uncorr(y) vs corr(x) [GeV/c]",
291  njptarray-1,jetptxbins,njptarray-1,jetptxbins);
292 
293 }//end function Set
294 
295 // ------------ method called once each job just before starting event loop ------------
296 // ------------ after setting histograms --------------------
297 // ------------ This function is needed to save histogram errors -----------------------
298 void
300 {
301  for (std::map<std::string,TH1D*>::const_iterator ih=histocontainer_.begin();
302  ih!= histocontainer_.end(); ++ih) {
303  TH1D *htemp = (TH1D*) ih->second;
304  htemp->Sumw2();
305  }
306 }//end function Sumw2
307 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
void Fill(edm::View< pat::Jet >::const_iterator &, std::string)
bool IsbTag(const pat::Jet &JetCand, const std::string &operpoint, const std::string &tagger) const
Definition: bJetSelector.cc:26
std::map< std::string, TH1D * > histocontainer_
size_type size() const
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
PatBTagCommonHistos(const edm::ParameterSet &)
std::string BTagdiscriminator_
std::map< std::string, TH2D * > h2_
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:85