ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/CmsHi/JetAnalysis/src/RecHitTreeProducer.cc
Revision: 1.1
Committed: Sat Oct 16 13:46:13 2010 UTC (14 years, 6 months ago) by yilmaz
Content type: text/plain
Branch: MAIN
Log Message:
Hit Tree producer

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     // $Id: RecHitTreeProducer.cc,v 1.1 2010/09/23 13:28:28 yilmaz Exp $
17     //
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     edm::Handle<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > > jets1;
108     edm::Handle<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > > jets2;
109    
110     // typedef edm::Handle<edm::SortedCollection<EcalRecHit,edm::StrictWeakOrdering<EcalRecHit> > >::const_iterator EcalIterator;
111     typedef vector<EcalRecHit>::const_iterator EcalIterator;
112    
113     edm::Handle<reco::CaloJetCollection> signalJets;
114    
115     MyRecHit hbheRecHit;
116     MyRecHit hfRecHit;
117    
118     TTree* hbheTree;
119     TTree* hfTree;
120    
121     double cone;
122    
123     edm::Service<TFileService> fs;
124     const CentralityBins * cbins_;
125     const CaloGeometry *geo;
126    
127     edm::InputTag HcalRecHitHFSrc_;
128     edm::InputTag HcalRecHitHBHESrc_;
129    
130    
131     };
132    
133     //
134     // constants, enums and typedefs
135     //
136    
137     //
138     // static data member definitions
139     //
140    
141     //
142     // constructors and destructor
143     //
144     RecHitTreeProducer::RecHitTreeProducer(const edm::ParameterSet& iConfig) :
145     cbins_(0),
146     geo(0),
147     cone(0.5)
148     {
149     //now do what ever initialization is needed
150     HcalRecHitHFSrc_ = iConfig.getUntrackedParameter<edm::InputTag>("hcalHFRecHitSrc",edm::InputTag("hfreco"));
151     HcalRecHitHBHESrc_ = iConfig.getUntrackedParameter<edm::InputTag>("hcalHBHERecHitSrc",edm::InputTag("hbhereco"));
152     }
153    
154    
155     RecHitTreeProducer::~RecHitTreeProducer()
156     {
157    
158     // do anything here that needs to be done at desctruction time
159     // (e.g. close files, deallocate resources etc.)
160    
161     }
162    
163    
164     //
165     // member functions
166     //
167    
168     // ------------ method called to for each event ------------
169     void
170     RecHitTreeProducer::analyze(const edm::Event& ev, const edm::EventSetup& iSetup)
171     {
172    
173     hfRecHit.n = 0;
174     hbheRecHit.n = 0;
175    
176     edm::Handle<HFRecHitCollection> hfHits;
177     edm::Handle<HBHERecHitCollection> hbheHits;
178    
179     ev.getByLabel(HcalRecHitHFSrc_,hfHits);
180     ev.getByLabel(HcalRecHitHBHESrc_,hbheHits);
181    
182    
183     if(0 && !cbins_) cbins_ = getCentralityBinsFromDB(iSetup);
184    
185     if(!geo){
186     edm::ESHandle<CaloGeometry> pGeo;
187     iSetup.get<CaloGeometryRecord>().get(pGeo);
188     geo = pGeo.product();
189     }
190    
191     for(unsigned int i = 0; i < hfHits->size(); ++i){
192     const HFRecHit & hit= (*hfHits)[i];
193     hfRecHit.e[hfRecHit.n] = hit.energy();
194     hfRecHit.et[hfRecHit.n] = getEt(hit.id(),hit.energy());
195     hfRecHit.eta[hfRecHit.n] = getEta(hit.id());
196     hfRecHit.phi[hfRecHit.n] = getPhi(hit.id());
197     // hfRecHit.ieta[hfRecHit.n] = hit.id().ieta();
198     // hfRecHit.iphi[hfRecHit.n] = hit.id().iphi();
199    
200     hfRecHit.n++;
201     }
202    
203     for(unsigned int i = 0; i < hbheHits->size(); ++i){
204     const HBHERecHit & hit= (*hbheHits)[i];
205     hbheRecHit.e[hfRecHit.n] = hit.energy();
206     hbheRecHit.et[hfRecHit.n] = getEt(hit.id(),hit.energy());
207     hbheRecHit.eta[hfRecHit.n] = getEta(hit.id());
208     hbheRecHit.phi[hfRecHit.n] = getPhi(hit.id());
209     // hbheRecHit.ieta[hfRecHit.n] = hit.id().ieta();
210     // hbheRecHit.iphi[hfRecHit.n] = hit.id().iphi();
211    
212     hbheRecHit.n++;
213     }
214    
215     hbheTree->Fill();
216     hfTree->Fill();
217    
218     }
219    
220    
221     // ------------ method called once each job just before starting event loop ------------
222     void
223     RecHitTreeProducer::beginJob()
224     {
225    
226     hbheTree = fs->make<TTree>("hbhe","");
227     hbheTree->Branch("n",&hbheRecHit.n,"n/I");
228     hbheTree->Branch("e",hbheRecHit.e,"e[n]/F");
229     hbheTree->Branch("et",hbheRecHit.et,"et[n]/F");
230     hbheTree->Branch("eta",hbheRecHit.eta,"eta[n]/F");
231     hbheTree->Branch("phi",hbheRecHit.phi,"phi[n]/F");
232    
233     hfTree = fs->make<TTree>("hf","");
234     hfTree->Branch("n",&hfRecHit.n,"n/I");
235     hfTree->Branch("e",hfRecHit.e,"e[n]/F");
236     hfTree->Branch("et",hfRecHit.et,"et[n]/F");
237     hfTree->Branch("eta",hfRecHit.eta,"eta[n]/F");
238     hfTree->Branch("phi",hfRecHit.phi,"phi[n]/F");
239    
240     }
241    
242     // ------------ method called once each job just after ending the event loop ------------
243     void
244     RecHitTreeProducer::endJob() {
245     }
246    
247     double RecHitTreeProducer::getEt(const DetId &id, double energy){
248     const GlobalPoint& pos=geo->getPosition(id);
249     double et = energy*sin(pos.theta());
250     return et;
251     }
252    
253     double RecHitTreeProducer::getEta(const DetId &id){
254     const GlobalPoint& pos=geo->getPosition(id);
255     double et = pos.eta();
256     return et;
257     }
258    
259     double RecHitTreeProducer::getPhi(const DetId &id){
260     const GlobalPoint& pos=geo->getPosition(id);
261     double et = pos.phi();
262     return et;
263     }
264    
265    
266    
267     //define this as a plug-in
268     DEFINE_FWK_MODULE(RecHitTreeProducer);