ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/CmsHi/JetAnalysis/src/RecHitTreeProducer.cc
Revision: 1.2
Committed: Mon Oct 18 16:13:37 2010 UTC (14 years, 6 months ago) by yilmaz
Content type: text/plain
Branch: MAIN
Changes since 1.1: +55 -15 lines
Log Message:
compilable analyzers

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     // Created: Tue Sep 7 11:38:19 EDT 2010
16 yilmaz 1.2 // $Id: RecHitTreeProducer.cc,v 1.1 2010/10/16 13:46:13 yilmaz Exp $
17 yilmaz 1.1 //
18     //
19    
20    
21     // system include files
22     #include <memory>
23     #include <vector>
24     #include <iostream>
25    
26     // user include files
27     #include "FWCore/Framework/interface/Frameworkfwd.h"
28     #include "FWCore/Framework/interface/EDAnalyzer.h"
29    
30     #include "FWCore/Framework/interface/Event.h"
31     #include "FWCore/Framework/interface/MakerMacros.h"
32     #include "FWCore/Framework/interface/ESHandle.h"
33     #include "FWCore/ParameterSet/interface/ParameterSet.h"
34    
35     #include "DataFormats/DetId/interface/DetId.h"
36     #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
37     #include "Geometry/Records/interface/IdealGeometryRecord.h"
38     #include "Geometry/Records/interface/CaloGeometryRecord.h"
39     #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
40    
41     #include "CommonTools/UtilAlgos/interface/TFileService.h"
42     #include "FWCore/ServiceRegistry/interface/Service.h"
43    
44     #include "DataFormats/Math/interface/deltaR.h"
45    
46     #include "DataFormats/Common/interface/Handle.h"
47     #include "DataFormats/FWLite/interface/ChainEvent.h"
48     #include "FWCore/Utilities/interface/InputTag.h"
49     #include "DataFormats/HeavyIonEvent/interface/CentralityBins.h"
50     #include "DataFormats/CaloTowers/interface/CaloTower.h"
51     #include "DataFormats/HeavyIonEvent/interface/Centrality.h"
52     #include "SimDataFormats/HiGenData/interface/GenHIEvent.h"
53     #include "DataFormats/PatCandidates/interface/Jet.h"
54     #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
55    
56     #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
57     #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
58     #include "DataFormats/HcalDetId/interface/HcalDetId.h"
59     #include "DataFormats/EcalDetId/interface/EcalDetIdCollections.h"
60     #include "DataFormats/DetId/interface/DetId.h"
61    
62    
63     #include "TNtuple.h"
64    
65     using namespace std;
66    
67     #define MAXHITS 100000
68    
69     struct MyRecHit{
70    
71     int n;
72    
73     float e[MAXHITS];
74     float et[MAXHITS];
75     float eta[MAXHITS];
76     float phi[MAXHITS];
77    
78     };
79    
80    
81    
82     //
83     // class declaration
84     //
85    
86     class RecHitTreeProducer : public edm::EDAnalyzer {
87     public:
88     explicit RecHitTreeProducer(const edm::ParameterSet&);
89     ~RecHitTreeProducer();
90     double getEt(const DetId &id, double energy);
91     double getEta(const DetId &id);
92     double getPhi(const DetId &id);
93    
94    
95     private:
96     virtual void beginJob() ;
97     virtual void analyze(const edm::Event&, const edm::EventSetup&);
98     virtual void endJob() ;
99    
100     // ----------member data ---------------------------
101    
102    
103     edm::Handle<reco::Centrality> cent;
104     edm::Handle<vector<double> > ktRhos;
105     edm::Handle<vector<double> > akRhos;
106    
107 yilmaz 1.2 edm::Handle<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > > ebHits;
108     edm::Handle<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > > eeHits;
109    
110     edm::Handle<HFRecHitCollection> hfHits;
111     edm::Handle<HBHERecHitCollection> hbheHits;
112 yilmaz 1.1
113     typedef vector<EcalRecHit>::const_iterator EcalIterator;
114    
115     edm::Handle<reco::CaloJetCollection> signalJets;
116    
117     MyRecHit hbheRecHit;
118     MyRecHit hfRecHit;
119 yilmaz 1.2 MyRecHit ebRecHit;
120     MyRecHit eeRecHit;
121 yilmaz 1.1
122     TTree* hbheTree;
123     TTree* hfTree;
124 yilmaz 1.2 TTree* ebTree;
125     TTree* eeTree;
126     double cone;
127 yilmaz 1.1
128     edm::Service<TFileService> fs;
129     const CentralityBins * cbins_;
130     const CaloGeometry *geo;
131    
132     edm::InputTag HcalRecHitHFSrc_;
133     edm::InputTag HcalRecHitHBHESrc_;
134 yilmaz 1.2 edm::InputTag EBSrc_;
135     edm::InputTag EESrc_;
136 yilmaz 1.1
137    
138     };
139    
140     //
141     // constants, enums and typedefs
142     //
143    
144     //
145     // static data member definitions
146     //
147    
148     //
149     // constructors and destructor
150     //
151     RecHitTreeProducer::RecHitTreeProducer(const edm::ParameterSet& iConfig) :
152     cbins_(0),
153     geo(0),
154     cone(0.5)
155     {
156     //now do what ever initialization is needed
157 yilmaz 1.2
158     EBSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("EBRecHitSrc",edm::InputTag("ecalRecHit","EcalRecHitsEB"));
159     EESrc_ = iConfig.getUntrackedParameter<edm::InputTag>("EERecHitSrc",edm::InputTag("ecalRecHit","EcalRecHitsEE"));
160 yilmaz 1.1 HcalRecHitHFSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("hcalHFRecHitSrc",edm::InputTag("hfreco"));
161     HcalRecHitHBHESrc_ = iConfig.getUntrackedParameter<edm::InputTag>("hcalHBHERecHitSrc",edm::InputTag("hbhereco"));
162     }
163    
164    
165     RecHitTreeProducer::~RecHitTreeProducer()
166     {
167    
168     // do anything here that needs to be done at desctruction time
169     // (e.g. close files, deallocate resources etc.)
170    
171     }
172    
173    
174     //
175     // member functions
176     //
177    
178     // ------------ method called to for each event ------------
179     void
180     RecHitTreeProducer::analyze(const edm::Event& ev, const edm::EventSetup& iSetup)
181     {
182    
183     hfRecHit.n = 0;
184     hbheRecHit.n = 0;
185 yilmaz 1.2 ebRecHit.n = 0;
186     eeRecHit.n = 0;
187 yilmaz 1.1
188 yilmaz 1.2 ev.getByLabel(EBSrc_,ebHits);
189     ev.getByLabel(EESrc_,eeHits);
190 yilmaz 1.1
191     ev.getByLabel(HcalRecHitHFSrc_,hfHits);
192     ev.getByLabel(HcalRecHitHBHESrc_,hbheHits);
193    
194     if(0 && !cbins_) cbins_ = getCentralityBinsFromDB(iSetup);
195    
196     if(!geo){
197     edm::ESHandle<CaloGeometry> pGeo;
198     iSetup.get<CaloGeometryRecord>().get(pGeo);
199     geo = pGeo.product();
200     }
201    
202     for(unsigned int i = 0; i < hfHits->size(); ++i){
203     const HFRecHit & hit= (*hfHits)[i];
204     hfRecHit.e[hfRecHit.n] = hit.energy();
205     hfRecHit.et[hfRecHit.n] = getEt(hit.id(),hit.energy());
206     hfRecHit.eta[hfRecHit.n] = getEta(hit.id());
207     hfRecHit.phi[hfRecHit.n] = getPhi(hit.id());
208     hfRecHit.n++;
209     }
210    
211     for(unsigned int i = 0; i < hbheHits->size(); ++i){
212     const HBHERecHit & hit= (*hbheHits)[i];
213     hbheRecHit.e[hfRecHit.n] = hit.energy();
214     hbheRecHit.et[hfRecHit.n] = getEt(hit.id(),hit.energy());
215     hbheRecHit.eta[hfRecHit.n] = getEta(hit.id());
216     hbheRecHit.phi[hfRecHit.n] = getPhi(hit.id());
217 yilmaz 1.2 hbheRecHit.n++;
218     }
219    
220     for(unsigned int i = 0; i < ebHits->size(); ++i){
221     const EcalRecHit & hit= (*ebHits)[i];
222     ebRecHit.e[ebRecHit.n] = hit.energy();
223     ebRecHit.et[ebRecHit.n] = getEt(hit.id(),hit.energy());
224     ebRecHit.eta[ebRecHit.n] = getEta(hit.id());
225     ebRecHit.phi[ebRecHit.n] = getPhi(hit.id());
226     ebRecHit.n++;
227     }
228 yilmaz 1.1
229 yilmaz 1.2 for(unsigned int i = 0; i < eeHits->size(); ++i){
230     const EcalRecHit & hit= (*eeHits)[i];
231     eeRecHit.e[eeRecHit.n] = hit.energy();
232     eeRecHit.et[eeRecHit.n] = getEt(hit.id(),hit.energy());
233     eeRecHit.eta[eeRecHit.n] = getEta(hit.id());
234     eeRecHit.phi[eeRecHit.n] = getPhi(hit.id());
235     eeRecHit.n++;
236 yilmaz 1.1 }
237    
238 yilmaz 1.2 eeTree->Fill();
239     ebTree->Fill();
240    
241 yilmaz 1.1 hbheTree->Fill();
242     hfTree->Fill();
243    
244     }
245    
246    
247     // ------------ method called once each job just before starting event loop ------------
248     void
249     RecHitTreeProducer::beginJob()
250     {
251    
252     hbheTree = fs->make<TTree>("hbhe","");
253     hbheTree->Branch("n",&hbheRecHit.n,"n/I");
254     hbheTree->Branch("e",hbheRecHit.e,"e[n]/F");
255     hbheTree->Branch("et",hbheRecHit.et,"et[n]/F");
256     hbheTree->Branch("eta",hbheRecHit.eta,"eta[n]/F");
257     hbheTree->Branch("phi",hbheRecHit.phi,"phi[n]/F");
258    
259     hfTree = fs->make<TTree>("hf","");
260     hfTree->Branch("n",&hfRecHit.n,"n/I");
261     hfTree->Branch("e",hfRecHit.e,"e[n]/F");
262     hfTree->Branch("et",hfRecHit.et,"et[n]/F");
263     hfTree->Branch("eta",hfRecHit.eta,"eta[n]/F");
264     hfTree->Branch("phi",hfRecHit.phi,"phi[n]/F");
265    
266 yilmaz 1.2 eeTree = fs->make<TTree>("ee","");
267     eeTree->Branch("n",&eeRecHit.n,"n/I");
268     eeTree->Branch("e",eeRecHit.e,"e[n]/F");
269     eeTree->Branch("et",eeRecHit.et,"et[n]/F");
270     eeTree->Branch("eta",eeRecHit.eta,"eta[n]/F");
271     eeTree->Branch("phi",eeRecHit.phi,"phi[n]/F");
272    
273     ebTree = fs->make<TTree>("eb","");
274     ebTree->Branch("n",&ebRecHit.n,"n/I");
275     ebTree->Branch("e",ebRecHit.e,"e[n]/F");
276     ebTree->Branch("et",ebRecHit.et,"et[n]/F");
277     ebTree->Branch("eta",ebRecHit.eta,"eta[n]/F");
278     ebTree->Branch("phi",ebRecHit.phi,"phi[n]/F");
279    
280 yilmaz 1.1 }
281    
282     // ------------ method called once each job just after ending the event loop ------------
283     void
284     RecHitTreeProducer::endJob() {
285     }
286    
287     double RecHitTreeProducer::getEt(const DetId &id, double energy){
288     const GlobalPoint& pos=geo->getPosition(id);
289     double et = energy*sin(pos.theta());
290     return et;
291     }
292    
293     double RecHitTreeProducer::getEta(const DetId &id){
294     const GlobalPoint& pos=geo->getPosition(id);
295     double et = pos.eta();
296     return et;
297     }
298    
299     double RecHitTreeProducer::getPhi(const DetId &id){
300     const GlobalPoint& pos=geo->getPosition(id);
301     double et = pos.phi();
302     return et;
303     }
304    
305    
306    
307     //define this as a plug-in
308     DEFINE_FWK_MODULE(RecHitTreeProducer);