ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/CmsHi/JetAnalysis/src/RecHitTreeProducer.cc
Revision: 1.17
Committed: Thu Nov 3 16:12:42 2011 UTC (13 years, 6 months ago) by yjlee
Content type: text/plain
Branch: MAIN
Changes since 1.16: +18 -11 lines
Log Message:
Add option to turn off HF hit

File Contents

# User Rev Content
1 yilmaz 1.1 // -*- C++ -*-
2     //
3     // Package: RecHitTreeProducer
4     // Class: RecHitTreeProducer
5     //
6     /**\class RecHitTreeProducer RecHitTreeProducer.cc CmsHi/RecHitTreeProducer/src/RecHitTreeProducer.cc
7    
8     Description: [one line class summary]
9    
10     Implementation:
11     [Notes on implementation]
12     */
13     //
14     // Original Author: Yetkin Yilmaz
15 yjlee 1.17 // Modified: Frank Ma, Yen-Jie Lee
16 yilmaz 1.1 // Created: Tue Sep 7 11:38:19 EDT 2010
17 yjlee 1.17 // $Id: RecHitTreeProducer.cc,v 1.16 2011/10/25 13:41:08 yjlee Exp $
18 yilmaz 1.1 //
19     //
20    
21 yjlee 1.16 #define versionTag "v1"
22 yilmaz 1.1 // system include files
23     #include <memory>
24     #include <vector>
25     #include <iostream>
26    
27     // user include files
28     #include "FWCore/Framework/interface/Frameworkfwd.h"
29     #include "FWCore/Framework/interface/EDAnalyzer.h"
30    
31     #include "FWCore/Framework/interface/Event.h"
32     #include "FWCore/Framework/interface/MakerMacros.h"
33     #include "FWCore/Framework/interface/ESHandle.h"
34     #include "FWCore/ParameterSet/interface/ParameterSet.h"
35    
36     #include "DataFormats/DetId/interface/DetId.h"
37     #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
38     #include "Geometry/Records/interface/IdealGeometryRecord.h"
39     #include "Geometry/Records/interface/CaloGeometryRecord.h"
40     #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
41    
42     #include "CommonTools/UtilAlgos/interface/TFileService.h"
43     #include "FWCore/ServiceRegistry/interface/Service.h"
44    
45     #include "DataFormats/Math/interface/deltaR.h"
46    
47     #include "DataFormats/Common/interface/Handle.h"
48     #include "DataFormats/FWLite/interface/ChainEvent.h"
49     #include "FWCore/Utilities/interface/InputTag.h"
50     #include "DataFormats/HeavyIonEvent/interface/CentralityBins.h"
51     #include "DataFormats/CaloTowers/interface/CaloTower.h"
52     #include "DataFormats/HeavyIonEvent/interface/Centrality.h"
53     #include "SimDataFormats/HiGenData/interface/GenHIEvent.h"
54     #include "DataFormats/PatCandidates/interface/Jet.h"
55     #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
56    
57     #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
58     #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
59     #include "DataFormats/HcalDetId/interface/HcalDetId.h"
60     #include "DataFormats/EcalDetId/interface/EcalDetIdCollections.h"
61 yilmaz 1.6 #include "DataFormats/EgammaReco/interface/BasicClusterFwd.h"
62 yilmaz 1.1 #include "DataFormats/DetId/interface/DetId.h"
63    
64 frankma 1.14 #include "DataFormats/VertexReco/interface/Vertex.h"
65     #include "DataFormats/VertexReco/interface/VertexFwd.h"
66    
67 yilmaz 1.1
68     #include "TNtuple.h"
69    
70     using namespace std;
71    
72 nart 1.3 #define MAXHITS 1000000
73 yilmaz 1.1
74     struct MyRecHit{
75 yilmaz 1.9 int depth[MAXHITS];
76 yilmaz 1.1 int n;
77    
78     float e[MAXHITS];
79     float et[MAXHITS];
80     float eta[MAXHITS];
81     float phi[MAXHITS];
82 nart 1.3 bool isjet[MAXHITS];
83 frankma 1.14 float etvtx[MAXHITS];
84     float etavtx[MAXHITS];
85     float emEtVtx[MAXHITS];
86     float hadEtVtx[MAXHITS];
87 yilmaz 1.1
88 yilmaz 1.6 float jtpt;
89     float jteta;
90     float jtphi;
91    
92 yilmaz 1.1 };
93    
94    
95 yilmaz 1.7 struct MyBkg{
96     int n;
97     float rho[50];
98     float sigma[50];
99     };
100    
101 yilmaz 1.1
102     //
103     // class declaration
104     //
105    
106     class RecHitTreeProducer : public edm::EDAnalyzer {
107     public:
108     explicit RecHitTreeProducer(const edm::ParameterSet&);
109     ~RecHitTreeProducer();
110     double getEt(const DetId &id, double energy);
111     double getEta(const DetId &id);
112     double getPhi(const DetId &id);
113 frankma 1.14 reco::Vertex::Point getVtx(const edm::Event& ev);
114    
115 yilmaz 1.1
116    
117     private:
118     virtual void beginJob() ;
119     virtual void analyze(const edm::Event&, const edm::EventSetup&);
120     virtual void endJob() ;
121    
122     // ----------member data ---------------------------
123    
124     edm::Handle<reco::Centrality> cent;
125     edm::Handle<vector<double> > ktRhos;
126     edm::Handle<vector<double> > akRhos;
127    
128 yilmaz 1.2 edm::Handle<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > > ebHits;
129     edm::Handle<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > > eeHits;
130    
131     edm::Handle<HFRecHitCollection> hfHits;
132     edm::Handle<HBHERecHitCollection> hbheHits;
133 yilmaz 1.1
134 yilmaz 1.6 edm::Handle<reco::BasicClusterCollection> bClusters;
135     edm::Handle<CaloTowerCollection> towers;
136 frankma 1.14 edm::Handle<reco::VertexCollection> vtxs;
137 yilmaz 1.6
138 nart 1.3 typedef vector<EcalRecHit>::const_iterator EcalIterator;
139    
140     edm::Handle<reco::CaloJetCollection> jets;
141 yilmaz 1.7
142     edm::Handle<std::vector<double> > rhos;
143     edm::Handle<std::vector<double> > sigmas;
144 nart 1.3
145 yilmaz 1.1 MyRecHit hbheRecHit;
146     MyRecHit hfRecHit;
147 yilmaz 1.2 MyRecHit ebRecHit;
148     MyRecHit eeRecHit;
149 yilmaz 1.6 MyRecHit myBC;
150     MyRecHit myTowers;
151 yilmaz 1.7 MyBkg bkg;
152    
153 yilmaz 1.9 TNtuple* nt;
154 yilmaz 1.1 TTree* hbheTree;
155     TTree* hfTree;
156 yilmaz 1.2 TTree* ebTree;
157     TTree* eeTree;
158 yilmaz 1.6 TTree* bcTree;
159     TTree* towerTree;
160 yilmaz 1.7 TTree* bkgTree;
161 yilmaz 1.6
162 yilmaz 1.2 double cone;
163 yilmaz 1.9 double hfTowerThreshold_;
164     double hfLongThreshold_;
165     double hfShortThreshold_;
166     double hbheThreshold_;
167     double ebThreshold_;
168     double eeThreshold_;
169 frankma 1.13
170     double hbhePtMin_;
171     double hfPtMin_;
172     double ebPtMin_;
173     double eePtMin_;
174     double towerPtMin_;
175    
176 yilmaz 1.1 edm::Service<TFileService> fs;
177     const CentralityBins * cbins_;
178     const CaloGeometry *geo;
179    
180     edm::InputTag HcalRecHitHFSrc_;
181     edm::InputTag HcalRecHitHBHESrc_;
182 yilmaz 1.2 edm::InputTag EBSrc_;
183     edm::InputTag EESrc_;
184 yilmaz 1.6 edm::InputTag BCSrc_;
185     edm::InputTag TowerSrc_;
186 frankma 1.14 edm::InputTag VtxSrc_;
187 yilmaz 1.6
188 nart 1.3 edm::InputTag JetSrc_;
189 yilmaz 1.6
190 yilmaz 1.7 edm::InputTag FastJetTag_;
191    
192 yilmaz 1.6 bool useJets_;
193     bool doBasicClusters_;
194 yilmaz 1.7 bool doTowers_;
195     bool doEcal_;
196     bool doHcal_;
197 yjlee 1.17 bool doHF_;
198 frankma 1.15 bool hasVtx_;
199 yilmaz 1.1
200 yilmaz 1.7 bool doFastJet_;
201 yilmaz 1.9
202     bool doEbyEonly_;
203 yilmaz 1.1 };
204    
205     //
206     // constants, enums and typedefs
207     //
208    
209     //
210     // static data member definitions
211     //
212    
213     //
214     // constructors and destructor
215     //
216     RecHitTreeProducer::RecHitTreeProducer(const edm::ParameterSet& iConfig) :
217     cbins_(0),
218     geo(0),
219     cone(0.5)
220     {
221     //now do what ever initialization is needed
222 yilmaz 1.2 EBSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("EBRecHitSrc",edm::InputTag("ecalRecHit","EcalRecHitsEB"));
223     EESrc_ = iConfig.getUntrackedParameter<edm::InputTag>("EERecHitSrc",edm::InputTag("ecalRecHit","EcalRecHitsEE"));
224 yilmaz 1.1 HcalRecHitHFSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("hcalHFRecHitSrc",edm::InputTag("hfreco"));
225     HcalRecHitHBHESrc_ = iConfig.getUntrackedParameter<edm::InputTag>("hcalHBHERecHitSrc",edm::InputTag("hbhereco"));
226 yilmaz 1.6 BCSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("BasicClusterSrc1",edm::InputTag("ecalRecHit","EcalRecHitsEB","RECO"));
227     TowerSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("towersSrc",edm::InputTag("towerMaker"));
228 frankma 1.14 VtxSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("vtxSrc",edm::InputTag("hiSelectedVertex"));
229 yilmaz 1.10 JetSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("JetSrc",edm::InputTag("iterativeConePu5CaloJets"));
230 yilmaz 1.7 useJets_ = iConfig.getUntrackedParameter<bool>("useJets",true);
231 yilmaz 1.6 doBasicClusters_ = iConfig.getUntrackedParameter<bool>("doBasicClusters",false);
232 yilmaz 1.7 doTowers_ = iConfig.getUntrackedParameter<bool>("doTowers",true);
233     doEcal_ = iConfig.getUntrackedParameter<bool>("doEcal",true);
234     doHcal_ = iConfig.getUntrackedParameter<bool>("doHcal",true);
235 yjlee 1.17 doHF_ = iConfig.getUntrackedParameter<bool>("doHF",true);
236 frankma 1.15 hasVtx_ = iConfig.getUntrackedParameter<bool>("hasVtx",false);
237 yilmaz 1.7 doFastJet_ = iConfig.getUntrackedParameter<bool>("doFastJet",true);
238     FastJetTag_ = iConfig.getUntrackedParameter<edm::InputTag>("FastJetTag",edm::InputTag("kt4CaloJets"));
239 yilmaz 1.9 doEbyEonly_ = iConfig.getUntrackedParameter<bool>("doEbyEonly",false);
240     hfTowerThreshold_ = iConfig.getUntrackedParameter<double>("HFtowerMin",3.);
241 yilmaz 1.11 hfLongThreshold_ = iConfig.getUntrackedParameter<double>("HFlongMin",0.5);
242     hfShortThreshold_ = iConfig.getUntrackedParameter<double>("HFshortMin",0.85);
243 frankma 1.13 hbhePtMin_ = iConfig.getUntrackedParameter<double>("HBHETreePtMin",0);
244     hfPtMin_ = iConfig.getUntrackedParameter<double>("HFTreePtMin",0);
245     ebPtMin_ = iConfig.getUntrackedParameter<double>("EBTreePtMin",0);
246     eePtMin_ = iConfig.getUntrackedParameter<double>("EETreePtMin",0.);
247     towerPtMin_ = iConfig.getUntrackedParameter<double>("TowerTreePtMin",0.);
248 yilmaz 1.1 }
249    
250    
251     RecHitTreeProducer::~RecHitTreeProducer()
252     {
253    
254     // do anything here that needs to be done at desctruction time
255     // (e.g. close files, deallocate resources etc.)
256    
257     }
258    
259    
260     //
261     // member functions
262     //
263    
264     // ------------ method called to for each event ------------
265     void
266     RecHitTreeProducer::analyze(const edm::Event& ev, const edm::EventSetup& iSetup)
267     {
268    
269     hfRecHit.n = 0;
270     hbheRecHit.n = 0;
271 yilmaz 1.2 ebRecHit.n = 0;
272     eeRecHit.n = 0;
273 yilmaz 1.7 myBC.n = 0;
274     myTowers.n = 0;
275     bkg.n = 0;
276 yilmaz 1.1
277 frankma 1.14 // get vertex
278 frankma 1.15 reco::Vertex::Point vtx;
279     if (hasVtx_) vtx = getVtx(ev);
280 frankma 1.14
281 yjlee 1.17 if(doEcal_){
282     ev.getByLabel(EBSrc_,ebHits);
283     ev.getByLabel(EESrc_,eeHits);
284 yilmaz 1.7 }
285 yjlee 1.17
286 yilmaz 1.7 if(doHcal_){
287 yjlee 1.17 ev.getByLabel(HcalRecHitHBHESrc_,hbheHits);
288     }
289     if(doHF_){
290     ev.getByLabel(HcalRecHitHFSrc_,hfHits);
291 yilmaz 1.7 }
292 yjlee 1.17
293 yilmaz 1.6 if(useJets_) {
294     ev.getByLabel(JetSrc_,jets);
295     }
296 nart 1.3
297 yilmaz 1.6 if(doBasicClusters_){
298     ev.getByLabel(BCSrc_,bClusters);
299 nart 1.3 }
300 yilmaz 1.6
301 yilmaz 1.7 if(doTowers_){
302     ev.getByLabel(TowerSrc_,towers);
303     }
304 yilmaz 1.6
305 yilmaz 1.7 if(doFastJet_){
306     ev.getByLabel(edm::InputTag(FastJetTag_.label(),"rhos",FastJetTag_.process()),rhos);
307     ev.getByLabel(edm::InputTag(FastJetTag_.label(),"sigmas",FastJetTag_.process()),sigmas);
308     bkg.n = rhos->size();
309     for(unsigned int i = 0; i < rhos->size(); ++i){
310     bkg.rho[i] = (*rhos)[i];
311     bkg.sigma[i] = (*sigmas)[i];
312     }
313     }
314 nart 1.3
315 yilmaz 1.6 if(0 && !cbins_) cbins_ = getCentralityBinsFromDB(iSetup);
316 yilmaz 1.1
317     if(!geo){
318     edm::ESHandle<CaloGeometry> pGeo;
319     iSetup.get<CaloGeometryRecord>().get(pGeo);
320     geo = pGeo.product();
321     }
322 yilmaz 1.9
323 yilmaz 1.12 int nHFlongPlus = 0;
324     int nHFshortPlus = 0;
325     int nHFtowerPlus = 0;
326     int nHFlongMinus = 0;
327     int nHFshortMinus = 0;
328     int nHFtowerMinus = 0;
329    
330    
331 nart 1.3
332 yjlee 1.17 if(doHF_){
333 yilmaz 1.1 for(unsigned int i = 0; i < hfHits->size(); ++i){
334     const HFRecHit & hit= (*hfHits)[i];
335     hfRecHit.e[hfRecHit.n] = hit.energy();
336     hfRecHit.et[hfRecHit.n] = getEt(hit.id(),hit.energy());
337     hfRecHit.eta[hfRecHit.n] = getEta(hit.id());
338     hfRecHit.phi[hfRecHit.n] = getPhi(hit.id());
339 nart 1.3 hfRecHit.isjet[hfRecHit.n] = false;
340 yilmaz 1.9 hfRecHit.depth[hfRecHit.n] = hit.id().depth();
341 yilmaz 1.12
342     if(hit.id().ieta() > 0){
343 yjlee 1.17 if(hit.energy() > hfShortThreshold_ && hit.id().depth() != 1) nHFshortPlus++;
344     if(hit.energy() > hfLongThreshold_ && hit.id().depth() == 1) nHFlongPlus++;
345 yilmaz 1.12 }else{
346     if(hit.energy() > hfShortThreshold_ && hit.id().depth() != 1) nHFshortMinus++;
347     if(hit.energy() > hfLongThreshold_ && hit.id().depth() == 1) nHFlongMinus++;
348     }
349 yilmaz 1.9
350 yilmaz 1.7 if(useJets_){
351 nart 1.3 for(unsigned int j = 0 ; j < jets->size(); ++j){
352     const reco::Jet& jet = (*jets)[j];
353     double dr = reco::deltaR(hfRecHit.eta[hfRecHit.n],hfRecHit.phi[hfRecHit.n],jet.eta(),jet.phi());
354     if(dr < cone){ hfRecHit.isjet[hfRecHit.n] = true; }
355     }
356     }
357 frankma 1.13 if (hfRecHit.et[hfRecHit.n]>=hfPtMin_) hfRecHit.n++;
358 yilmaz 1.1 }
359 yjlee 1.17
360     if(doHcal_ && !doEbyEonly_){
361 yilmaz 1.1 for(unsigned int i = 0; i < hbheHits->size(); ++i){
362     const HBHERecHit & hit= (*hbheHits)[i];
363 frankma 1.13 if (getEt(hit.id(),hit.energy())<hbhePtMin_) continue;
364 nart 1.3 hbheRecHit.e[hbheRecHit.n] = hit.energy();
365     hbheRecHit.et[hbheRecHit.n] = getEt(hit.id(),hit.energy());
366     hbheRecHit.eta[hbheRecHit.n] = getEta(hit.id());
367     hbheRecHit.phi[hbheRecHit.n] = getPhi(hit.id());
368     hbheRecHit.isjet[hbheRecHit.n] = false;
369 yilmaz 1.7 if(useJets_){
370 nart 1.3 for(unsigned int j = 0 ; j < jets->size(); ++j){
371     const reco::Jet& jet = (*jets)[j];
372     double dr = reco::deltaR(hbheRecHit.eta[hbheRecHit.n],hbheRecHit.phi[hbheRecHit.n],jet.eta(),jet.phi());
373     if(dr < cone){ hbheRecHit.isjet[hbheRecHit.n] = true; }
374     }
375     }
376 yilmaz 1.2 hbheRecHit.n++;
377     }
378 yilmaz 1.7 }
379 yilmaz 1.9 }
380     if(doEcal_ && !doEbyEonly_){
381 yilmaz 1.2 for(unsigned int i = 0; i < ebHits->size(); ++i){
382     const EcalRecHit & hit= (*ebHits)[i];
383 frankma 1.13 if (getEt(hit.id(),hit.energy())<ebPtMin_) continue;
384 yilmaz 1.2 ebRecHit.e[ebRecHit.n] = hit.energy();
385     ebRecHit.et[ebRecHit.n] = getEt(hit.id(),hit.energy());
386     ebRecHit.eta[ebRecHit.n] = getEta(hit.id());
387     ebRecHit.phi[ebRecHit.n] = getPhi(hit.id());
388 nart 1.3 ebRecHit.isjet[ebRecHit.n] = false;
389 yilmaz 1.7 if(useJets_){
390 nart 1.3 for(unsigned int j = 0 ; j < jets->size(); ++j){
391     const reco::Jet& jet = (*jets)[j];
392     double dr = reco::deltaR(ebRecHit.eta[ebRecHit.n],ebRecHit.phi[ebRecHit.n],jet.eta(),jet.phi());
393     if(dr < cone){ ebRecHit.isjet[ebRecHit.n] = true; }
394     }
395     }
396 yilmaz 1.2 ebRecHit.n++;
397     }
398 nart 1.3
399 yilmaz 1.2 for(unsigned int i = 0; i < eeHits->size(); ++i){
400     const EcalRecHit & hit= (*eeHits)[i];
401 frankma 1.13 if (getEt(hit.id(),hit.energy())<eePtMin_) continue;
402 yilmaz 1.2 eeRecHit.e[eeRecHit.n] = hit.energy();
403     eeRecHit.et[eeRecHit.n] = getEt(hit.id(),hit.energy());
404     eeRecHit.eta[eeRecHit.n] = getEta(hit.id());
405     eeRecHit.phi[eeRecHit.n] = getPhi(hit.id());
406 nart 1.3 eeRecHit.isjet[eeRecHit.n] = false;
407 yilmaz 1.7 if(useJets_){
408 nart 1.3 for(unsigned int j = 0 ; j < jets->size(); ++j){
409     const reco::Jet& jet = (*jets)[j];
410     double dr = reco::deltaR(eeRecHit.eta[eeRecHit.n],eeRecHit.phi[eeRecHit.n],jet.eta(),jet.phi());
411     if(dr < cone){ eeRecHit.isjet[eeRecHit.n] = true; }
412     }
413     }
414 yilmaz 1.2 eeRecHit.n++;
415 yilmaz 1.1 }
416 yilmaz 1.7 }
417    
418     if(doTowers_){
419    
420     for(unsigned int i = 0; i < towers->size(); ++i){
421     const CaloTower & hit= (*towers)[i];
422 frankma 1.13 if (getEt(hit.id(),hit.energy())<towerPtMin_) continue;
423 yilmaz 1.7 myTowers.e[myTowers.n] = hit.energy();
424     myTowers.et[myTowers.n] = getEt(hit.id(),hit.energy());
425     myTowers.eta[myTowers.n] = getEta(hit.id());
426     myTowers.phi[myTowers.n] = getPhi(hit.id());
427     myTowers.isjet[myTowers.n] = false;
428 frankma 1.15 if (hasVtx_) {
429     myTowers.etvtx[myTowers.n] = hit.p4(vtx).Et();
430     myTowers.etavtx[myTowers.n] = hit.p4(vtx).Eta();
431     myTowers.emEtVtx[myTowers.n] = hit.emEt(vtx);
432     myTowers.hadEtVtx[myTowers.n] = hit.hadEt(vtx);
433     }
434 yilmaz 1.9
435 yilmaz 1.12 if(hit.ieta() > 29 && hit.energy() > hfTowerThreshold_) nHFtowerPlus++;
436     if(hit.ieta() < -29 && hit.energy() > hfTowerThreshold_) nHFtowerMinus++;
437 yilmaz 1.9
438 yilmaz 1.7 if(useJets_){
439     for(unsigned int j = 0 ; j < jets->size(); ++j){
440     const reco::Jet& jet = (*jets)[j];
441     double dr = reco::deltaR(myTowers.eta[myTowers.n],myTowers.phi[myTowers.n],jet.eta(),jet.phi());
442     if(dr < cone){ myTowers.isjet[myTowers.n] = true; }
443     }
444     }
445     myTowers.n++;
446     }
447 yilmaz 1.6
448 yilmaz 1.7 }
449 yilmaz 1.6
450 yilmaz 1.9 if(doBasicClusters_ && !doEbyEonly_){
451 yilmaz 1.6 for(unsigned int j = 0 ; j < jets->size(); ++j){
452     const reco::Jet& jet = (*jets)[j];
453     myBC.n = 0;
454     myBC.jtpt = jet.pt();
455     myBC.jteta = jet.eta();
456     myBC.jtphi = jet.phi();
457    
458     for(unsigned int i = 0; i < bClusters->size(); ++i){
459     const reco::BasicCluster & bc= (*bClusters)[i];
460     double dr = reco::deltaR(bc.eta(),bc.phi(),jet.eta(),jet.phi());
461     if(dr < cone){
462     myBC.e[myBC.n] = bc.energy();
463     myBC.et[myBC.n] = bc.energy()*sin(bc.position().theta());
464     myBC.eta[myBC.n] = bc.eta();
465     myBC.phi[myBC.n] = bc.phi();
466     myBC.n++;
467     }
468     }
469     bcTree->Fill();
470     }
471     }
472 yilmaz 1.7
473 yilmaz 1.9 if(!doEbyEonly_){
474     towerTree->Fill();
475    
476     eeTree->Fill();
477     ebTree->Fill();
478    
479     hbheTree->Fill();
480     hfTree->Fill();
481    
482     if (doFastJet_) {
483     bkgTree->Fill();
484     }
485     }
486    
487 yilmaz 1.12 nt->Fill(nHFtowerPlus,nHFtowerMinus,nHFlongPlus,nHFlongMinus,nHFshortPlus,nHFshortMinus);
488 nart 1.3
489 yilmaz 1.1 }
490    
491    
492     // ------------ method called once each job just before starting event loop ------------
493     void
494     RecHitTreeProducer::beginJob()
495     {
496 nart 1.3
497 yjlee 1.16 hbheTree = fs->make<TTree>("hbhe",versionTag);
498 yilmaz 1.1 hbheTree->Branch("n",&hbheRecHit.n,"n/I");
499     hbheTree->Branch("e",hbheRecHit.e,"e[n]/F");
500     hbheTree->Branch("et",hbheRecHit.et,"et[n]/F");
501     hbheTree->Branch("eta",hbheRecHit.eta,"eta[n]/F");
502     hbheTree->Branch("phi",hbheRecHit.phi,"phi[n]/F");
503 yilmaz 1.5 hbheTree->Branch("isjet",hbheRecHit.isjet,"isjet[n]/O");
504 nart 1.3
505 yjlee 1.16 hfTree = fs->make<TTree>("hf",versionTag);
506 yilmaz 1.1 hfTree->Branch("n",&hfRecHit.n,"n/I");
507     hfTree->Branch("e",hfRecHit.e,"e[n]/F");
508     hfTree->Branch("et",hfRecHit.et,"et[n]/F");
509     hfTree->Branch("eta",hfRecHit.eta,"eta[n]/F");
510     hfTree->Branch("phi",hfRecHit.phi,"phi[n]/F");
511 yilmaz 1.9 hfTree->Branch("depth",hfRecHit.depth,"depth[n]/I");
512 yilmaz 1.5 hfTree->Branch("isjet",hfRecHit.isjet,"isjet[n]/O");
513 yilmaz 1.1
514 yjlee 1.16 eeTree = fs->make<TTree>("ee",versionTag);
515 yilmaz 1.2 eeTree->Branch("n",&eeRecHit.n,"n/I");
516     eeTree->Branch("e",eeRecHit.e,"e[n]/F");
517     eeTree->Branch("et",eeRecHit.et,"et[n]/F");
518     eeTree->Branch("eta",eeRecHit.eta,"eta[n]/F");
519     eeTree->Branch("phi",eeRecHit.phi,"phi[n]/F");
520 yilmaz 1.5 eeTree->Branch("isjet",eeRecHit.isjet,"isjet[n]/O");
521 nart 1.3
522 yjlee 1.16 ebTree = fs->make<TTree>("eb",versionTag);
523 yilmaz 1.2 ebTree->Branch("n",&ebRecHit.n,"n/I");
524     ebTree->Branch("e",ebRecHit.e,"e[n]/F");
525     ebTree->Branch("et",ebRecHit.et,"et[n]/F");
526     ebTree->Branch("eta",ebRecHit.eta,"eta[n]/F");
527     ebTree->Branch("phi",ebRecHit.phi,"phi[n]/F");
528 yilmaz 1.4 ebTree->Branch("isjet",ebRecHit.isjet,"isjet[n]/O");
529 yilmaz 1.2
530 yjlee 1.16 towerTree = fs->make<TTree>("tower",versionTag);
531 yilmaz 1.7 towerTree->Branch("n",&myTowers.n,"n/I");
532     towerTree->Branch("e",myTowers.e,"e[n]/F");
533     towerTree->Branch("et",myTowers.et,"et[n]/F");
534     towerTree->Branch("eta",myTowers.eta,"eta[n]/F");
535     towerTree->Branch("phi",myTowers.phi,"phi[n]/F");
536     towerTree->Branch("isjet",myTowers.isjet,"isjet[n]/O");
537 frankma 1.15 if (hasVtx_) {
538     towerTree->Branch("etvtx",myTowers.etvtx,"etvtx[n]/F");
539     towerTree->Branch("etavtx",myTowers.etavtx,"etavtx[n]/F");
540     towerTree->Branch("emEtVtx",myTowers.emEtVtx,"emEtVtx[n]/F");
541     towerTree->Branch("hadEtVtx",myTowers.hadEtVtx,"hadEtVtx[n]/F");
542     }
543 yilmaz 1.7
544    
545 yilmaz 1.6 if(doBasicClusters_){
546 yjlee 1.16 bcTree = fs->make<TTree>("bc",versionTag);
547 yilmaz 1.6 bcTree->Branch("n",&myBC.n,"n/I");
548     bcTree->Branch("e",myBC.e,"e[n]/F");
549     bcTree->Branch("et",myBC.et,"et[n]/F");
550     bcTree->Branch("eta",myBC.eta,"eta[n]/F");
551     bcTree->Branch("phi",myBC.phi,"phi[n]/F");
552     bcTree->Branch("jtpt",&myBC.jtpt,"jtpt/F");
553     bcTree->Branch("jteta",&myBC.jteta,"jteta/F");
554     bcTree->Branch("jtphi",&myBC.jtphi,"jtphi/F");
555     // bcTree->Branch("isjet",bcRecHit.isjet,"isjet[n]/O");
556     }
557    
558 yilmaz 1.7 if(doFastJet_){
559 yjlee 1.16 bkgTree = fs->make<TTree>("bkg",versionTag);
560 yilmaz 1.7 bkgTree->Branch("n",&bkg.n,"n/I");
561     bkgTree->Branch("rho",bkg.rho,"rho[n]/F");
562     bkgTree->Branch("sigma",bkg.sigma,"sigma[n]/F");
563     }
564 yilmaz 1.9
565 yilmaz 1.12 nt = fs->make<TNtuple>("ntEvent","","nHFplus:nHFminus:nHFlongPlus:nHFlongMinus:nHFshortPlus:nHFshortMinus");
566 yilmaz 1.6
567 yilmaz 1.1 }
568    
569     // ------------ method called once each job just after ending the event loop ------------
570     void
571     RecHitTreeProducer::endJob() {
572     }
573    
574     double RecHitTreeProducer::getEt(const DetId &id, double energy){
575     const GlobalPoint& pos=geo->getPosition(id);
576     double et = energy*sin(pos.theta());
577     return et;
578     }
579 nart 1.3
580 yilmaz 1.1 double RecHitTreeProducer::getEta(const DetId &id){
581     const GlobalPoint& pos=geo->getPosition(id);
582     double et = pos.eta();
583     return et;
584     }
585    
586     double RecHitTreeProducer::getPhi(const DetId &id){
587     const GlobalPoint& pos=geo->getPosition(id);
588     double et = pos.phi();
589     return et;
590     }
591    
592 frankma 1.14 reco::Vertex::Point RecHitTreeProducer::getVtx(const edm::Event& ev)
593     {
594     ev.getByLabel(VtxSrc_,vtxs);
595     int greatestvtx = 0;
596     int nVertex = vtxs->size();
597    
598     for (unsigned int i = 0 ; i< vtxs->size(); ++i){
599     unsigned int daughter = (*vtxs)[i].tracksSize();
600     if( daughter > (*vtxs)[greatestvtx].tracksSize()) greatestvtx = i;
601     //cout <<"Vertex: "<< (*vtxs)[i].position().z()<<" "<<daughter<<endl;
602     }
603    
604     if(nVertex<=0){
605     return reco::Vertex::Point(-999,-999,-999);
606     }
607     return (*vtxs)[greatestvtx].position();
608     }
609 yilmaz 1.1
610     //define this as a plug-in
611     DEFINE_FWK_MODULE(RecHitTreeProducer);