ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/CmsHi/JetAnalysis/src/RecHitTreeProducer.cc
Revision: 1.16
Committed: Tue Oct 25 13:41:08 2011 UTC (13 years, 6 months ago) by yjlee
Content type: text/plain
Branch: MAIN
Changes since 1.15: +9 -9 lines
Log Message:
Add version tag

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 frankma 1.14 // Modified: Frank Ma
16 yilmaz 1.1 // Created: Tue Sep 7 11:38:19 EDT 2010
17 yjlee 1.16 // $Id: RecHitTreeProducer.cc,v 1.15 2011/10/14 15:06:14 frankma 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 frankma 1.15 bool hasVtx_;
198 yilmaz 1.1
199 yilmaz 1.7 bool doFastJet_;
200 yilmaz 1.9
201     bool doEbyEonly_;
202 yilmaz 1.1 };
203    
204     //
205     // constants, enums and typedefs
206     //
207    
208     //
209     // static data member definitions
210     //
211    
212     //
213     // constructors and destructor
214     //
215     RecHitTreeProducer::RecHitTreeProducer(const edm::ParameterSet& iConfig) :
216     cbins_(0),
217     geo(0),
218     cone(0.5)
219     {
220     //now do what ever initialization is needed
221 yilmaz 1.2 EBSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("EBRecHitSrc",edm::InputTag("ecalRecHit","EcalRecHitsEB"));
222     EESrc_ = iConfig.getUntrackedParameter<edm::InputTag>("EERecHitSrc",edm::InputTag("ecalRecHit","EcalRecHitsEE"));
223 yilmaz 1.1 HcalRecHitHFSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("hcalHFRecHitSrc",edm::InputTag("hfreco"));
224     HcalRecHitHBHESrc_ = iConfig.getUntrackedParameter<edm::InputTag>("hcalHBHERecHitSrc",edm::InputTag("hbhereco"));
225 yilmaz 1.6 BCSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("BasicClusterSrc1",edm::InputTag("ecalRecHit","EcalRecHitsEB","RECO"));
226     TowerSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("towersSrc",edm::InputTag("towerMaker"));
227 frankma 1.14 VtxSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("vtxSrc",edm::InputTag("hiSelectedVertex"));
228 yilmaz 1.10 JetSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("JetSrc",edm::InputTag("iterativeConePu5CaloJets"));
229 yilmaz 1.7 useJets_ = iConfig.getUntrackedParameter<bool>("useJets",true);
230 yilmaz 1.6 doBasicClusters_ = iConfig.getUntrackedParameter<bool>("doBasicClusters",false);
231 yilmaz 1.7 doTowers_ = iConfig.getUntrackedParameter<bool>("doTowers",true);
232     doEcal_ = iConfig.getUntrackedParameter<bool>("doEcal",true);
233     doHcal_ = iConfig.getUntrackedParameter<bool>("doHcal",true);
234 frankma 1.15 hasVtx_ = iConfig.getUntrackedParameter<bool>("hasVtx",false);
235 yilmaz 1.7 doFastJet_ = iConfig.getUntrackedParameter<bool>("doFastJet",true);
236     FastJetTag_ = iConfig.getUntrackedParameter<edm::InputTag>("FastJetTag",edm::InputTag("kt4CaloJets"));
237 yilmaz 1.9 doEbyEonly_ = iConfig.getUntrackedParameter<bool>("doEbyEonly",false);
238     hfTowerThreshold_ = iConfig.getUntrackedParameter<double>("HFtowerMin",3.);
239 yilmaz 1.11 hfLongThreshold_ = iConfig.getUntrackedParameter<double>("HFlongMin",0.5);
240     hfShortThreshold_ = iConfig.getUntrackedParameter<double>("HFshortMin",0.85);
241 frankma 1.13 hbhePtMin_ = iConfig.getUntrackedParameter<double>("HBHETreePtMin",0);
242     hfPtMin_ = iConfig.getUntrackedParameter<double>("HFTreePtMin",0);
243     ebPtMin_ = iConfig.getUntrackedParameter<double>("EBTreePtMin",0);
244     eePtMin_ = iConfig.getUntrackedParameter<double>("EETreePtMin",0.);
245     towerPtMin_ = iConfig.getUntrackedParameter<double>("TowerTreePtMin",0.);
246 yilmaz 1.1 }
247    
248    
249     RecHitTreeProducer::~RecHitTreeProducer()
250     {
251    
252     // do anything here that needs to be done at desctruction time
253     // (e.g. close files, deallocate resources etc.)
254    
255     }
256    
257    
258     //
259     // member functions
260     //
261    
262     // ------------ method called to for each event ------------
263     void
264     RecHitTreeProducer::analyze(const edm::Event& ev, const edm::EventSetup& iSetup)
265     {
266    
267     hfRecHit.n = 0;
268     hbheRecHit.n = 0;
269 yilmaz 1.2 ebRecHit.n = 0;
270     eeRecHit.n = 0;
271 yilmaz 1.7 myBC.n = 0;
272     myTowers.n = 0;
273     bkg.n = 0;
274 yilmaz 1.1
275 frankma 1.14 // get vertex
276 frankma 1.15 reco::Vertex::Point vtx;
277     if (hasVtx_) vtx = getVtx(ev);
278 frankma 1.14
279 yilmaz 1.7 if(doEcal_){
280 yilmaz 1.2 ev.getByLabel(EBSrc_,ebHits);
281     ev.getByLabel(EESrc_,eeHits);
282 yilmaz 1.7 }
283     if(doHcal_){
284 yilmaz 1.1 ev.getByLabel(HcalRecHitHFSrc_,hfHits);
285     ev.getByLabel(HcalRecHitHBHESrc_,hbheHits);
286 yilmaz 1.7 }
287 yilmaz 1.6 if(useJets_) {
288     ev.getByLabel(JetSrc_,jets);
289     }
290 nart 1.3
291 yilmaz 1.6 if(doBasicClusters_){
292     ev.getByLabel(BCSrc_,bClusters);
293 nart 1.3 }
294 yilmaz 1.6
295 yilmaz 1.7 if(doTowers_){
296     ev.getByLabel(TowerSrc_,towers);
297     }
298 yilmaz 1.6
299 yilmaz 1.7 if(doFastJet_){
300     ev.getByLabel(edm::InputTag(FastJetTag_.label(),"rhos",FastJetTag_.process()),rhos);
301     ev.getByLabel(edm::InputTag(FastJetTag_.label(),"sigmas",FastJetTag_.process()),sigmas);
302     bkg.n = rhos->size();
303     for(unsigned int i = 0; i < rhos->size(); ++i){
304     bkg.rho[i] = (*rhos)[i];
305     bkg.sigma[i] = (*sigmas)[i];
306     }
307     }
308 nart 1.3
309 yilmaz 1.6 if(0 && !cbins_) cbins_ = getCentralityBinsFromDB(iSetup);
310 yilmaz 1.1
311     if(!geo){
312     edm::ESHandle<CaloGeometry> pGeo;
313     iSetup.get<CaloGeometryRecord>().get(pGeo);
314     geo = pGeo.product();
315     }
316 yilmaz 1.9
317 yilmaz 1.12 int nHFlongPlus = 0;
318     int nHFshortPlus = 0;
319     int nHFtowerPlus = 0;
320     int nHFlongMinus = 0;
321     int nHFshortMinus = 0;
322     int nHFtowerMinus = 0;
323    
324    
325 nart 1.3
326 yilmaz 1.7 if(doHcal_){
327 yilmaz 1.1 for(unsigned int i = 0; i < hfHits->size(); ++i){
328     const HFRecHit & hit= (*hfHits)[i];
329     hfRecHit.e[hfRecHit.n] = hit.energy();
330     hfRecHit.et[hfRecHit.n] = getEt(hit.id(),hit.energy());
331     hfRecHit.eta[hfRecHit.n] = getEta(hit.id());
332     hfRecHit.phi[hfRecHit.n] = getPhi(hit.id());
333 nart 1.3 hfRecHit.isjet[hfRecHit.n] = false;
334 yilmaz 1.9 hfRecHit.depth[hfRecHit.n] = hit.id().depth();
335 yilmaz 1.12
336     if(hit.id().ieta() > 0){
337     if(hit.energy() > hfShortThreshold_ && hit.id().depth() != 1) nHFshortPlus++;
338     if(hit.energy() > hfLongThreshold_ && hit.id().depth() == 1) nHFlongPlus++;
339     }else{
340     if(hit.energy() > hfShortThreshold_ && hit.id().depth() != 1) nHFshortMinus++;
341     if(hit.energy() > hfLongThreshold_ && hit.id().depth() == 1) nHFlongMinus++;
342     }
343 yilmaz 1.9
344 yilmaz 1.7 if(useJets_){
345 nart 1.3 for(unsigned int j = 0 ; j < jets->size(); ++j){
346     const reco::Jet& jet = (*jets)[j];
347     double dr = reco::deltaR(hfRecHit.eta[hfRecHit.n],hfRecHit.phi[hfRecHit.n],jet.eta(),jet.phi());
348     if(dr < cone){ hfRecHit.isjet[hfRecHit.n] = true; }
349     }
350     }
351 frankma 1.13 if (hfRecHit.et[hfRecHit.n]>=hfPtMin_) hfRecHit.n++;
352 yilmaz 1.1 }
353 yilmaz 1.9 if(!doEbyEonly_){
354 yilmaz 1.1 for(unsigned int i = 0; i < hbheHits->size(); ++i){
355     const HBHERecHit & hit= (*hbheHits)[i];
356 frankma 1.13 if (getEt(hit.id(),hit.energy())<hbhePtMin_) continue;
357 nart 1.3 hbheRecHit.e[hbheRecHit.n] = hit.energy();
358     hbheRecHit.et[hbheRecHit.n] = getEt(hit.id(),hit.energy());
359     hbheRecHit.eta[hbheRecHit.n] = getEta(hit.id());
360     hbheRecHit.phi[hbheRecHit.n] = getPhi(hit.id());
361     hbheRecHit.isjet[hbheRecHit.n] = false;
362 yilmaz 1.7 if(useJets_){
363 nart 1.3 for(unsigned int j = 0 ; j < jets->size(); ++j){
364     const reco::Jet& jet = (*jets)[j];
365     double dr = reco::deltaR(hbheRecHit.eta[hbheRecHit.n],hbheRecHit.phi[hbheRecHit.n],jet.eta(),jet.phi());
366     if(dr < cone){ hbheRecHit.isjet[hbheRecHit.n] = true; }
367     }
368     }
369 yilmaz 1.2 hbheRecHit.n++;
370     }
371 yilmaz 1.7 }
372 yilmaz 1.9 }
373     if(doEcal_ && !doEbyEonly_){
374 yilmaz 1.2 for(unsigned int i = 0; i < ebHits->size(); ++i){
375     const EcalRecHit & hit= (*ebHits)[i];
376 frankma 1.13 if (getEt(hit.id(),hit.energy())<ebPtMin_) continue;
377 yilmaz 1.2 ebRecHit.e[ebRecHit.n] = hit.energy();
378     ebRecHit.et[ebRecHit.n] = getEt(hit.id(),hit.energy());
379     ebRecHit.eta[ebRecHit.n] = getEta(hit.id());
380     ebRecHit.phi[ebRecHit.n] = getPhi(hit.id());
381 nart 1.3 ebRecHit.isjet[ebRecHit.n] = false;
382 yilmaz 1.7 if(useJets_){
383 nart 1.3 for(unsigned int j = 0 ; j < jets->size(); ++j){
384     const reco::Jet& jet = (*jets)[j];
385     double dr = reco::deltaR(ebRecHit.eta[ebRecHit.n],ebRecHit.phi[ebRecHit.n],jet.eta(),jet.phi());
386     if(dr < cone){ ebRecHit.isjet[ebRecHit.n] = true; }
387     }
388     }
389 yilmaz 1.2 ebRecHit.n++;
390     }
391 nart 1.3
392 yilmaz 1.2 for(unsigned int i = 0; i < eeHits->size(); ++i){
393     const EcalRecHit & hit= (*eeHits)[i];
394 frankma 1.13 if (getEt(hit.id(),hit.energy())<eePtMin_) continue;
395 yilmaz 1.2 eeRecHit.e[eeRecHit.n] = hit.energy();
396     eeRecHit.et[eeRecHit.n] = getEt(hit.id(),hit.energy());
397     eeRecHit.eta[eeRecHit.n] = getEta(hit.id());
398     eeRecHit.phi[eeRecHit.n] = getPhi(hit.id());
399 nart 1.3 eeRecHit.isjet[eeRecHit.n] = false;
400 yilmaz 1.7 if(useJets_){
401 nart 1.3 for(unsigned int j = 0 ; j < jets->size(); ++j){
402     const reco::Jet& jet = (*jets)[j];
403     double dr = reco::deltaR(eeRecHit.eta[eeRecHit.n],eeRecHit.phi[eeRecHit.n],jet.eta(),jet.phi());
404     if(dr < cone){ eeRecHit.isjet[eeRecHit.n] = true; }
405     }
406     }
407 yilmaz 1.2 eeRecHit.n++;
408 yilmaz 1.1 }
409 yilmaz 1.7 }
410    
411     if(doTowers_){
412    
413     for(unsigned int i = 0; i < towers->size(); ++i){
414     const CaloTower & hit= (*towers)[i];
415 frankma 1.13 if (getEt(hit.id(),hit.energy())<towerPtMin_) continue;
416 yilmaz 1.7 myTowers.e[myTowers.n] = hit.energy();
417     myTowers.et[myTowers.n] = getEt(hit.id(),hit.energy());
418     myTowers.eta[myTowers.n] = getEta(hit.id());
419     myTowers.phi[myTowers.n] = getPhi(hit.id());
420     myTowers.isjet[myTowers.n] = false;
421 frankma 1.15 if (hasVtx_) {
422     myTowers.etvtx[myTowers.n] = hit.p4(vtx).Et();
423     myTowers.etavtx[myTowers.n] = hit.p4(vtx).Eta();
424     myTowers.emEtVtx[myTowers.n] = hit.emEt(vtx);
425     myTowers.hadEtVtx[myTowers.n] = hit.hadEt(vtx);
426     }
427 yilmaz 1.9
428 yilmaz 1.12 if(hit.ieta() > 29 && hit.energy() > hfTowerThreshold_) nHFtowerPlus++;
429     if(hit.ieta() < -29 && hit.energy() > hfTowerThreshold_) nHFtowerMinus++;
430 yilmaz 1.9
431 yilmaz 1.7 if(useJets_){
432     for(unsigned int j = 0 ; j < jets->size(); ++j){
433     const reco::Jet& jet = (*jets)[j];
434     double dr = reco::deltaR(myTowers.eta[myTowers.n],myTowers.phi[myTowers.n],jet.eta(),jet.phi());
435     if(dr < cone){ myTowers.isjet[myTowers.n] = true; }
436     }
437     }
438     myTowers.n++;
439     }
440 yilmaz 1.6
441 yilmaz 1.7 }
442 yilmaz 1.6
443 yilmaz 1.9 if(doBasicClusters_ && !doEbyEonly_){
444 yilmaz 1.6 for(unsigned int j = 0 ; j < jets->size(); ++j){
445     const reco::Jet& jet = (*jets)[j];
446     myBC.n = 0;
447     myBC.jtpt = jet.pt();
448     myBC.jteta = jet.eta();
449     myBC.jtphi = jet.phi();
450    
451     for(unsigned int i = 0; i < bClusters->size(); ++i){
452     const reco::BasicCluster & bc= (*bClusters)[i];
453     double dr = reco::deltaR(bc.eta(),bc.phi(),jet.eta(),jet.phi());
454     if(dr < cone){
455     myBC.e[myBC.n] = bc.energy();
456     myBC.et[myBC.n] = bc.energy()*sin(bc.position().theta());
457     myBC.eta[myBC.n] = bc.eta();
458     myBC.phi[myBC.n] = bc.phi();
459     myBC.n++;
460     }
461     }
462     bcTree->Fill();
463     }
464     }
465 yilmaz 1.7
466 yilmaz 1.9 if(!doEbyEonly_){
467     towerTree->Fill();
468    
469     eeTree->Fill();
470     ebTree->Fill();
471    
472     hbheTree->Fill();
473     hfTree->Fill();
474    
475     if (doFastJet_) {
476     bkgTree->Fill();
477     }
478     }
479    
480 yilmaz 1.12 nt->Fill(nHFtowerPlus,nHFtowerMinus,nHFlongPlus,nHFlongMinus,nHFshortPlus,nHFshortMinus);
481 nart 1.3
482 yilmaz 1.1 }
483    
484    
485     // ------------ method called once each job just before starting event loop ------------
486     void
487     RecHitTreeProducer::beginJob()
488     {
489 nart 1.3
490 yjlee 1.16 hbheTree = fs->make<TTree>("hbhe",versionTag);
491 yilmaz 1.1 hbheTree->Branch("n",&hbheRecHit.n,"n/I");
492     hbheTree->Branch("e",hbheRecHit.e,"e[n]/F");
493     hbheTree->Branch("et",hbheRecHit.et,"et[n]/F");
494     hbheTree->Branch("eta",hbheRecHit.eta,"eta[n]/F");
495     hbheTree->Branch("phi",hbheRecHit.phi,"phi[n]/F");
496 yilmaz 1.5 hbheTree->Branch("isjet",hbheRecHit.isjet,"isjet[n]/O");
497 nart 1.3
498 yjlee 1.16 hfTree = fs->make<TTree>("hf",versionTag);
499 yilmaz 1.1 hfTree->Branch("n",&hfRecHit.n,"n/I");
500     hfTree->Branch("e",hfRecHit.e,"e[n]/F");
501     hfTree->Branch("et",hfRecHit.et,"et[n]/F");
502     hfTree->Branch("eta",hfRecHit.eta,"eta[n]/F");
503     hfTree->Branch("phi",hfRecHit.phi,"phi[n]/F");
504 yilmaz 1.9 hfTree->Branch("depth",hfRecHit.depth,"depth[n]/I");
505 yilmaz 1.5 hfTree->Branch("isjet",hfRecHit.isjet,"isjet[n]/O");
506 yilmaz 1.1
507 yjlee 1.16 eeTree = fs->make<TTree>("ee",versionTag);
508 yilmaz 1.2 eeTree->Branch("n",&eeRecHit.n,"n/I");
509     eeTree->Branch("e",eeRecHit.e,"e[n]/F");
510     eeTree->Branch("et",eeRecHit.et,"et[n]/F");
511     eeTree->Branch("eta",eeRecHit.eta,"eta[n]/F");
512     eeTree->Branch("phi",eeRecHit.phi,"phi[n]/F");
513 yilmaz 1.5 eeTree->Branch("isjet",eeRecHit.isjet,"isjet[n]/O");
514 nart 1.3
515 yjlee 1.16 ebTree = fs->make<TTree>("eb",versionTag);
516 yilmaz 1.2 ebTree->Branch("n",&ebRecHit.n,"n/I");
517     ebTree->Branch("e",ebRecHit.e,"e[n]/F");
518     ebTree->Branch("et",ebRecHit.et,"et[n]/F");
519     ebTree->Branch("eta",ebRecHit.eta,"eta[n]/F");
520     ebTree->Branch("phi",ebRecHit.phi,"phi[n]/F");
521 yilmaz 1.4 ebTree->Branch("isjet",ebRecHit.isjet,"isjet[n]/O");
522 yilmaz 1.2
523 yjlee 1.16 towerTree = fs->make<TTree>("tower",versionTag);
524 yilmaz 1.7 towerTree->Branch("n",&myTowers.n,"n/I");
525     towerTree->Branch("e",myTowers.e,"e[n]/F");
526     towerTree->Branch("et",myTowers.et,"et[n]/F");
527     towerTree->Branch("eta",myTowers.eta,"eta[n]/F");
528     towerTree->Branch("phi",myTowers.phi,"phi[n]/F");
529     towerTree->Branch("isjet",myTowers.isjet,"isjet[n]/O");
530 frankma 1.15 if (hasVtx_) {
531     towerTree->Branch("etvtx",myTowers.etvtx,"etvtx[n]/F");
532     towerTree->Branch("etavtx",myTowers.etavtx,"etavtx[n]/F");
533     towerTree->Branch("emEtVtx",myTowers.emEtVtx,"emEtVtx[n]/F");
534     towerTree->Branch("hadEtVtx",myTowers.hadEtVtx,"hadEtVtx[n]/F");
535     }
536 yilmaz 1.7
537    
538 yilmaz 1.6 if(doBasicClusters_){
539 yjlee 1.16 bcTree = fs->make<TTree>("bc",versionTag);
540 yilmaz 1.6 bcTree->Branch("n",&myBC.n,"n/I");
541     bcTree->Branch("e",myBC.e,"e[n]/F");
542     bcTree->Branch("et",myBC.et,"et[n]/F");
543     bcTree->Branch("eta",myBC.eta,"eta[n]/F");
544     bcTree->Branch("phi",myBC.phi,"phi[n]/F");
545     bcTree->Branch("jtpt",&myBC.jtpt,"jtpt/F");
546     bcTree->Branch("jteta",&myBC.jteta,"jteta/F");
547     bcTree->Branch("jtphi",&myBC.jtphi,"jtphi/F");
548     // bcTree->Branch("isjet",bcRecHit.isjet,"isjet[n]/O");
549     }
550    
551 yilmaz 1.7 if(doFastJet_){
552 yjlee 1.16 bkgTree = fs->make<TTree>("bkg",versionTag);
553 yilmaz 1.7 bkgTree->Branch("n",&bkg.n,"n/I");
554     bkgTree->Branch("rho",bkg.rho,"rho[n]/F");
555     bkgTree->Branch("sigma",bkg.sigma,"sigma[n]/F");
556     }
557 yilmaz 1.9
558 yilmaz 1.12 nt = fs->make<TNtuple>("ntEvent","","nHFplus:nHFminus:nHFlongPlus:nHFlongMinus:nHFshortPlus:nHFshortMinus");
559 yilmaz 1.6
560 yilmaz 1.1 }
561    
562     // ------------ method called once each job just after ending the event loop ------------
563     void
564     RecHitTreeProducer::endJob() {
565     }
566    
567     double RecHitTreeProducer::getEt(const DetId &id, double energy){
568     const GlobalPoint& pos=geo->getPosition(id);
569     double et = energy*sin(pos.theta());
570     return et;
571     }
572 nart 1.3
573 yilmaz 1.1 double RecHitTreeProducer::getEta(const DetId &id){
574     const GlobalPoint& pos=geo->getPosition(id);
575     double et = pos.eta();
576     return et;
577     }
578    
579     double RecHitTreeProducer::getPhi(const DetId &id){
580     const GlobalPoint& pos=geo->getPosition(id);
581     double et = pos.phi();
582     return et;
583     }
584    
585 frankma 1.14 reco::Vertex::Point RecHitTreeProducer::getVtx(const edm::Event& ev)
586     {
587     ev.getByLabel(VtxSrc_,vtxs);
588     int greatestvtx = 0;
589     int nVertex = vtxs->size();
590    
591     for (unsigned int i = 0 ; i< vtxs->size(); ++i){
592     unsigned int daughter = (*vtxs)[i].tracksSize();
593     if( daughter > (*vtxs)[greatestvtx].tracksSize()) greatestvtx = i;
594     //cout <<"Vertex: "<< (*vtxs)[i].position().z()<<" "<<daughter<<endl;
595     }
596    
597     if(nVertex<=0){
598     return reco::Vertex::Point(-999,-999,-999);
599     }
600     return (*vtxs)[greatestvtx].position();
601     }
602 yilmaz 1.1
603     //define this as a plug-in
604     DEFINE_FWK_MODULE(RecHitTreeProducer);