CMS 3D CMS Logo

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