CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
SUSYDQMAnalyzer.cc
Go to the documentation of this file.
1 //authors: Francesco Costanza (DESY)
2 // Dirk Kruecker (DESY)
3 //date: 05/05/11
4 
5 //================================================================
6 // CMS FW and Data Handlers
7 
9 
11 
18 
19 //================================================================
20 // Jet & Jet collections // MET & MET collections
21 
26 
31 //#include "DataFormats/JetReco/interface/JPTJet.h"
32 //#include "DataFormats/JetReco/interface/JPTJetCollection.h"
33 
41 
42 //================================================================
43 // SUSY Classes
44 
49 
50 //================================================================
51 // ROOT Classes
52 
53 #include "TH1.h"
54 #include "TVector2.h"
55 #include "TLorentzVector.h"
56 
57 //================================================================
58 // Ordinary C++ stuff
59 
60 #include <cmath>
61 #include <vector>
62 #include <fstream>
63 #include <sstream>
64 #include <string>
65 
66 using namespace edm;
67 using namespace reco;
68 using namespace math;
69 using namespace std;
70 
72  iConfig = pSet;
73 
74  SUSYFolder = iConfig.getParameter<std::string>("folderName");
75  // Load parameters
76  thePFMETCollectionToken =
77  consumes<reco::PFMETCollection>(iConfig.getParameter<edm::InputTag>("PFMETCollectionLabel"));
78  theCaloMETCollectionToken =
79  consumes<reco::CaloMETCollection>(iConfig.getParameter<edm::InputTag>("CaloMETCollectionLabel"));
80 
81  //remove TCMET and JPT related variables due to anticipated changes in RECO
82  //theTCMETCollectionToken = consumes<reco::METCollection> (iConfig.getParameter<edm::InputTag>("TCMETCollectionLabel"));
83 
84  theCaloJetCollectionToken =
85  consumes<reco::CaloJetCollection>(iConfig.getParameter<edm::InputTag>("CaloJetCollectionLabel"));
86  //theJPTJetCollectionToken = consumes<reco::JPTJetCollection> (iConfig.getParameter<edm::InputTag>("JPTJetCollectionLabel"));
87  thePFJetCollectionToken =
88  consumes<std::vector<reco::PFJet> >(iConfig.getParameter<edm::InputTag>("PFJetCollectionLabel"));
89 
90  _ptThreshold = iConfig.getParameter<double>("ptThreshold");
91  _maxNJets = iConfig.getParameter<double>("maxNJets");
92  _maxAbsEta = iConfig.getParameter<double>("maxAbsEta");
93 }
94 
95 const char* SUSYDQMAnalyzer::messageLoggerCatregory = "SUSYDQM";
96 
98  // if( dqm ) {
99  //===========================================================
100  // book HT histos.
101  std::string dir = SUSYFolder;
102  dir += "HT";
103  ibooker.setCurrentFolder(dir);
104  hCaloHT = ibooker.book1D("Calo_HT", "", 500, 0., 2000);
105  hPFHT = ibooker.book1D("PF_HT", "", 500, 0., 2000);
106  //hJPTHT = ibooker.book1D("JPT_HT" , "", 500, 0., 2000);
107  //===========================================================
108  // book MET histos.
109 
110  dir = SUSYFolder;
111  dir += "MET";
112  ibooker.setCurrentFolder(dir);
113  hCaloMET = ibooker.book1D("Calo_MET", "", 500, 0., 1000);
114  hPFMET = ibooker.book1D("PF_MET", "", 500, 0., 1000);
115  //hTCMET = ibooker.book1D("TC_MET" , "", 500, 0., 1000);
116 
117  //===========================================================
118  // book MHT histos.
119 
120  dir = SUSYFolder;
121  dir += "MHT";
122  ibooker.setCurrentFolder(dir);
123  hCaloMHT = ibooker.book1D("Calo_MHT", "", 500, 0., 1000);
124  hPFMHT = ibooker.book1D("PF_MHT", "", 500, 0., 1000);
125  //hJPTMHT = ibooker.book1D("JPT_MHT" , "", 500, 0., 1000);
126 
127  //===========================================================
128  // book alpha_T histos.
129 
130  dir = SUSYFolder;
131  dir += "Alpha_T";
132  ibooker.setCurrentFolder(dir);
133  hCaloAlpha_T = ibooker.book1D("Calo_AlphaT", "", 100, 0., 1.);
134  //hJPTAlpha_T = ibooker.book1D("PF_AlphaT" , "", 100, 0., 1.);
135  hPFAlpha_T = ibooker.book1D("PF_AlphaT", "", 100, 0., 1.);
136  // }
137 }
138 
140  //###########################################################
141  // HTand MHT
142 
143  //===========================================================
144  // Calo HT, MHT and alpha_T
145 
147 
148  iEvent.getByToken(theCaloJetCollectionToken, CaloJetcoll);
149 
150  if (!CaloJetcoll.isValid())
151  return;
152 
153  std::vector<math::XYZTLorentzVector> Ps;
154  for (reco::CaloJetCollection::const_iterator jet = CaloJetcoll->begin(); jet != CaloJetcoll->end(); ++jet) {
155  if ((jet->pt() > _ptThreshold) && (abs(jet->eta()) < _maxAbsEta)) {
156  if (Ps.size() > _maxNJets) {
157  edm::LogInfo(messageLoggerCatregory) << "NMax Jets exceded..";
158  break;
159  }
160  Ps.push_back(jet->p4());
161  }
162  }
163 
164  hCaloAlpha_T->Fill(alpha_T()(Ps));
165 
166  HT<reco::CaloJetCollection> CaloHT(CaloJetcoll, _ptThreshold, _maxAbsEta);
167 
168  hCaloHT->Fill(CaloHT.ScalarSum);
169  hCaloMHT->Fill(CaloHT.v.Mod());
170 
171  //===========================================================
172  // PF HT, MHT and alpha_T
173 
175 
176  iEvent.getByToken(thePFJetCollectionToken, PFjetcoll);
177 
178  if (!PFjetcoll.isValid())
179  return;
180 
181  Ps.clear();
182  for (reco::PFJetCollection::const_iterator jet = PFjetcoll->begin(); jet != PFjetcoll->end(); ++jet) {
183  if ((jet->pt() > _ptThreshold) && (abs(jet->eta()) < _maxAbsEta)) {
184  if (Ps.size() > _maxNJets) {
185  edm::LogInfo(messageLoggerCatregory) << "NMax Jets exceded..";
186  break;
187  }
188  Ps.push_back(jet->p4());
189  }
190  }
191  hPFAlpha_T->Fill(alpha_T()(Ps));
192 
193  HT<reco::PFJetCollection> PFHT(PFjetcoll, _ptThreshold, _maxAbsEta);
194 
195  hPFHT->Fill(PFHT.ScalarSum);
196  hPFMHT->Fill(PFHT.v.Mod());
197 
198  //===========================================================
199  // JPT HT, MHT and alpha_T
200 
201  //edm::Handle<reco::JPTJetCollection> JPTjetcoll;
202 
203  //iEvent.getByToken(theJPTJetCollectionToken, JPTjetcoll);
204 
205  //if(!JPTjetcoll.isValid()) return;
206 
207  //Ps.clear();
208  //for (reco::JPTJetCollection::const_iterator jet = JPTjetcoll->begin(); jet!=JPTjetcoll->end(); ++jet){
209  //if ((jet->pt()>_ptThreshold) && (abs(jet->eta())<_maxAbsEta)){
210  // if(Ps.size()>_maxNJets) {
211  // edm::LogInfo(messageLoggerCatregory)<<"NMax Jets exceded..";
212  // break;
213  // }
214  // Ps.push_back(jet->p4());
215  //}
216  //}
217  //hJPTAlpha_T->Fill( alpha_T()(Ps));
218 
219  //HT<reco::JPTJetCollection> JPTHT(JPTjetcoll, _ptThreshold, _maxAbsEta);
220 
221  //hJPTHT->Fill(JPTHT.ScalarSum);
222  //hJPTMHT->Fill(JPTHT.v.Mod());
223 
224  //###########################################################
225  // MET
226 
227  //===========================================================
228  // Calo MET
229 
231  iEvent.getByToken(theCaloMETCollectionToken, calometcoll);
232 
233  if (!calometcoll.isValid())
234  return;
235 
236  const CaloMETCollection* calometcol = calometcoll.product();
237  const CaloMET* calomet;
238  calomet = &(calometcol->front());
239 
240  hCaloMET->Fill(calomet->pt());
241 
242  //===========================================================
243  // PF MET
244 
246  iEvent.getByToken(thePFMETCollectionToken, pfmetcoll);
247 
248  if (!pfmetcoll.isValid())
249  return;
250 
251  const PFMETCollection* pfmetcol = pfmetcoll.product();
252  const PFMET* pfmet;
253  pfmet = &(pfmetcol->front());
254 
255  hPFMET->Fill(pfmet->pt());
256 
257  //===========================================================
258  // TC MET
259 
260  //edm::Handle<reco::METCollection> tcmetcoll;
261  //iEvent.getByToken(theTCMETCollectionToken, tcmetcoll);
262 
263  //if(!tcmetcoll.isValid()) return;
264 
265  //const METCollection *tcmetcol = tcmetcoll.product();
266  //const MET *tcmet;
267  //tcmet = &(tcmetcol->front());
268 
269  //hTCMET->Fill(tcmet->pt());
270 }
271 
~SUSYDQMAnalyzer() override
double pt() const final
transverse momentum
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
double ScalarSum
Definition: HT.h:28
Collection of Calo MET.
int iEvent
Definition: GenABIO.cc:224
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool isValid() const
Definition: HandleBase.h:70
Log< level::Info, false > LogInfo
T const * product() const
Definition: Handle.h:70
static const char * messageLoggerCatregory
SUSYDQMAnalyzer(const edm::ParameterSet &)
void analyze(const edm::Event &, const edm::EventSetup &) override
TVector2 v
Definition: HT.h:27
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
Definition: HT.h:21
Definition: Run.h:45
Collection of PF MET.