CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CentralityProducer.cc
Go to the documentation of this file.
1 //
2 // Original Author: Yetkin Yilmaz, Young Soo Park
3 // Created: Wed Jun 11 15:31:41 CEST 2008
4 //
5 //
6 
7 
8 // system include files
9 #include <memory>
10 #include <iostream>
11 
12 // user include files
17 
22 
24 
29 
38 
44 
45 using namespace std;
46 using namespace edm;
47 using namespace reco;
48 
49 //
50 // class declaration
51 //
52 
53 namespace reco {
55  public:
56  explicit CentralityProducer(const edm::ParameterSet&);
58 
59  private:
60  virtual void beginJob() override ;
61  virtual void produce(edm::Event&, const edm::EventSetup&) override;
62  virtual void endJob() override ;
63 
64  // ----------member data ---------------------------
65 
66  bool recoLevel_;
67 
75  bool reuseAny_;
77 
79 
81  double trackPtCut_;
82  double trackEtaCut_;
83  double hfEtaCut_;
84 
86 
97 
100 
103 
104 };
105 
106 //
107 // constants, enums and typedefs
108 //
109 
110 
111 //
112 // static data member definitions
113 //
114 
115 //
116 // constructors and destructor
117 //
118  CentralityProducer::CentralityProducer(const edm::ParameterSet& iConfig) {
119  //register your products
120  produceHFhits_ = iConfig.getParameter<bool>("produceHFhits");
121  produceHFtowers_ = iConfig.getParameter<bool>("produceHFtowers");
122  produceEcalhits_ = iConfig.getParameter<bool>("produceEcalhits");
123  produceZDChits_ = iConfig.getParameter<bool>("produceZDChits");
124  produceETmidRap_ = iConfig.getParameter<bool>("produceETmidRapidity");
125  producePixelhits_ = iConfig.getParameter<bool>("producePixelhits");
126  produceTracks_ = iConfig.getParameter<bool>("produceTracks");
127  producePixelTracks_ = iConfig.getParameter<bool>("producePixelTracks");
128 
129  midRapidityRange_ = iConfig.getParameter<double>("midRapidityRange");
130  trackPtCut_ = iConfig.getParameter<double>("trackPtCut");
131  trackEtaCut_ = iConfig.getParameter<double>("trackEtaCut");
132 
133  hfEtaCut_ = iConfig.getParameter<double>("hfEtaCut");
134 
135  if(produceHFhits_) srcHFhits_ = consumes<HFRecHitCollection>(iConfig.getParameter<edm::InputTag>("srcHFhits"));
136  if(produceHFtowers_ || produceETmidRap_) srcTowers_ = consumes<CaloTowerCollection>(iConfig.getParameter<edm::InputTag>("srcTowers"));
137 
138  if(produceEcalhits_){
139  srcEBhits_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("srcEBhits"));
140  srcEEhits_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("srcEEhits"));
141  }
142  if(produceZDChits_){
143  srcZDChits_ = consumes<ZDCRecHitCollection>(iConfig.getParameter<edm::InputTag>("srcZDChits"));
144  lowGainZDC_ = iConfig.getParameter<bool>("lowGainZDC");
145  }
146  if(producePixelhits_){
147  srcPixelhits_ = consumes<SiPixelRecHitCollection>(iConfig.getParameter<edm::InputTag>("srcPixelhits"));
148  doPixelCut_ = iConfig.getParameter<bool>("doPixelCut");
149  srcVertex_ = consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("srcVertex"));
150  }
151  if(produceTracks_) {
152  srcTracks_ = consumes<TrackCollection>(iConfig.getParameter<edm::InputTag>("srcTracks"));
153  srcVertex_ = consumes<VertexCollection>(iConfig.getParameter<edm::InputTag>("srcVertex"));
154  }
155  if(producePixelTracks_) srcPixelTracks_ = consumes<TrackCollection>(iConfig.getParameter<edm::InputTag>("srcPixelTracks"));
156 
157  reuseAny_ = iConfig.getParameter<bool>("reUseCentrality");
158  if(reuseAny_) reuseTag_ = consumes<Centrality>(iConfig.getParameter<edm::InputTag>("srcReUse"));
159 
160  useQuality_ = iConfig.getParameter<bool>("UseQuality");
161  trackQuality_ = TrackBase::qualityByName(iConfig.getParameter<std::string>("TrackQuality"));
162 
163  produces<Centrality>();
164 
165 }
166 
167 
168 CentralityProducer::~CentralityProducer()
169 {
170 
171  // do anything here that needs to be done at desctruction time
172  // (e.g. close files, deallocate resources etc.)
173 
174 }
175 
176 
177 //
178 // member functions
179 //
180 
181 // ------------ method called to produce the data ------------
182 void
183 CentralityProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
184 {
185 
186  if(producePixelhits_) iSetup.get<TrackerDigiGeometryRecord>().get(tGeo);
187  if(produceEcalhits_) iSetup.get<CaloGeometryRecord>().get(cGeo);
188 
189  std::auto_ptr<Centrality> creco(new Centrality());
190  Handle<Centrality> inputCentrality;
191 
192  if(reuseAny_) iEvent.getByToken(reuseTag_,inputCentrality);
193 
194  if(produceHFhits_){
195  creco->etHFhitSumPlus_ = 0;
196  creco->etHFhitSumMinus_ = 0;
197 
199  iEvent.getByToken(srcHFhits_,hits);
200  for( size_t ihit = 0; ihit<hits->size(); ++ ihit){
201  const HFRecHit & rechit = (*hits)[ ihit ];
202  if(rechit.id().ieta() > 0 )
203  creco->etHFhitSumPlus_ += rechit.energy();
204  if(rechit.id().ieta() < 0)
205  creco->etHFhitSumMinus_ += rechit.energy();
206  }
207  }else{
208  if(reuseAny_){
209  creco->etHFhitSumMinus_ = inputCentrality->EtHFhitSumMinus();
210  creco->etHFhitSumPlus_ = inputCentrality->EtHFhitSumPlus();
211  }
212  }
213 
214  if(produceHFtowers_ || produceETmidRap_){
215  creco->etHFtowerSumPlus_ = 0;
216  creco->etHFtowerSumMinus_ = 0;
217  creco->etMidRapiditySum_ = 0;
218 
220 
221  iEvent.getByToken(srcTowers_,towers);
222 
223  for( size_t i = 0; i<towers->size(); ++ i){
224  const CaloTower & tower = (*towers)[ i ];
225  double eta = tower.eta();
226  bool isHF = tower.ietaAbs() > 29;
227  if(produceHFtowers_){
228  if(isHF && eta > 0){
229  creco->etHFtowerSumPlus_ += tower.pt();
230  if(eta > hfEtaCut_) creco->etHFtruncatedPlus_ += tower.pt();
231  }
232  if(isHF && eta < 0){
233  creco->etHFtowerSumMinus_ += tower.pt();
234  if(eta < -hfEtaCut_) creco->etHFtruncatedMinus_ += tower.pt();
235  }
236  }else{
237  if(reuseAny_){
238  creco->etHFtowerSumMinus_ = inputCentrality->EtHFtowerSumMinus();
239  creco->etHFtowerSumPlus_ = inputCentrality->EtHFtowerSumPlus();
240  creco->etHFtruncatedMinus_ = inputCentrality->EtHFtruncatedMinus();
241  creco->etHFtruncatedPlus_ = inputCentrality->EtHFtruncatedPlus();
242  }
243  }
244  if(produceETmidRap_){
245  if(fabs(eta) < midRapidityRange_) creco->etMidRapiditySum_ += tower.pt()/(midRapidityRange_*2.);
246  }else if(reuseAny_) creco->etMidRapiditySum_ = inputCentrality->EtMidRapiditySum();
247  }
248  }else{
249  if(reuseAny_){
250  creco->etHFtowerSumMinus_ = inputCentrality->EtHFtowerSumMinus();
251  creco->etHFtowerSumPlus_ = inputCentrality->EtHFtowerSumPlus();
252  creco->etMidRapiditySum_ = inputCentrality->EtMidRapiditySum();
253  }
254  }
255 
256  if(produceEcalhits_){
257  creco->etEESumPlus_ = 0;
258  creco->etEESumMinus_ = 0;
259  creco->etEBSum_ = 0;
260 
263 
264  iEvent.getByToken(srcEBhits_,ebHits);
265  iEvent.getByToken(srcEEhits_,eeHits);
266 
267  for(unsigned int i = 0; i < ebHits->size(); ++i){
268  const EcalRecHit & hit= (*ebHits)[i];
269  const GlobalPoint& pos=cGeo->getPosition(hit.id());
270  double et = hit.energy()*sin(pos.theta());
271  creco->etEBSum_ += et;
272  }
273 
274  for(unsigned int i = 0; i < eeHits->size(); ++i){
275  const EcalRecHit & hit= (*eeHits)[i];
276  const GlobalPoint& pos=cGeo->getPosition(hit.id());
277  double et = hit.energy()*sin(pos.theta());
278  if(pos.z() > 0){
279  creco->etEESumPlus_ += et;
280  }else{
281  creco->etEESumMinus_ += et;
282  }
283  }
284  }else{
285  if(reuseAny_){
286  creco->etEESumMinus_ = inputCentrality->EtEESumMinus();
287  creco->etEESumPlus_ = inputCentrality->EtEESumPlus();
288  creco->etEBSum_ = inputCentrality->EtEBSum();
289  }
290  }
291 
292  if(producePixelhits_){
293  creco->pixelMultiplicity_ = 0;
296  iEvent.getByToken(srcPixelhits_,rchts);
297  rechits = rchts.product();
298  int nPixel =0 ;
299  for (SiPixelRecHitCollection::const_iterator it = rechits->begin(); it!=rechits->end();it++)
300  {
302  DetId detId = DetId(hits.detId());
303  SiPixelRecHitCollection::const_iterator recHitMatch = rechits->find(detId);
304  const SiPixelRecHitCollection::DetSet recHitRange = *recHitMatch;
305  for ( SiPixelRecHitCollection::DetSet::const_iterator recHitIterator = recHitRange.begin();
306  recHitIterator != recHitRange.end(); ++recHitIterator) {
307  // add selection if needed, now all hits.
308  if(doPixelCut_){
309  const SiPixelRecHit * recHit = &(*recHitIterator);
310  const PixelGeomDetUnit* pixelLayer = dynamic_cast<const PixelGeomDetUnit*> (tGeo->idToDet(recHit->geographicalId()));
311  GlobalPoint gpos = pixelLayer->toGlobal(recHit->localPosition());
312  math::XYZVector rechitPos(gpos.x(),gpos.y(),gpos.z());
313  double abeta = fabs(rechitPos.eta());
314  int clusterSize = recHit->cluster()->size();
315  if ( abeta < 0.5 && clusterSize < 1) continue;
316  if ( abeta > 0.5 && abeta < 1 && clusterSize < 2) continue;
317  if ( abeta > 1. && abeta < 1.5 && clusterSize < 3) continue;
318  if ( abeta > 1.5 && abeta < 2. && clusterSize < 4) continue;
319  if ( abeta > 2. && abeta < 2.5 && clusterSize < 6) continue;
320  if ( abeta > 2.5 && abeta < 5 && clusterSize < 9) continue;
321  }
322  nPixel++;
323  }
324  }
325  creco->pixelMultiplicity_ = nPixel;
326  }else{
327  if(reuseAny_){
328  creco->pixelMultiplicity_ = inputCentrality->multiplicityPixel();
329  }
330  }
331 
332  if(produceTracks_){
333 
334  double vx=-999.;
335  double vy=-999.;
336  double vz=-999.;
337  double vxError=-999.;
338  double vyError=-999.;
339  double vzError=-999.;
340  math::XYZVector vtxPos(0,0,0);
341 
342  Handle<VertexCollection> recoVertices;
343  iEvent.getByToken(srcVertex_,recoVertices);
344  unsigned int daughter = 0;
345  int greatestvtx = 0;
346 
347  for (unsigned int i = 0 ; i< recoVertices->size(); ++i){
348  daughter = (*recoVertices)[i].tracksSize();
349  if( daughter > (*recoVertices)[greatestvtx].tracksSize()) greatestvtx = i;
350  }
351 
352  if(recoVertices->size()>0){
353  vx = (*recoVertices)[greatestvtx].position().x();
354  vy = (*recoVertices)[greatestvtx].position().y();
355  vz = (*recoVertices)[greatestvtx].position().z();
356  vxError = (*recoVertices)[greatestvtx].xError();
357  vyError = (*recoVertices)[greatestvtx].yError();
358  vzError = (*recoVertices)[greatestvtx].zError();
359  }
360 
361  vtxPos = math::XYZVector(vx,vy,vz);
362 
364  iEvent.getByToken(srcTracks_,tracks);
365  int nTracks = 0;
366 
367  double trackCounter = 0;
368  double trackCounterEta = 0;
369  double trackCounterEtaPt = 0;
370 
371  for(unsigned int i = 0 ; i < tracks->size(); ++i){
372  const Track& track = (*tracks)[i];
373  if(useQuality_ && !track.quality(trackQuality_)) continue;
374 
375  if( track.pt() > trackPtCut_) trackCounter++;
376  if(fabs(track.eta())<trackEtaCut_) {
377  trackCounterEta++;
378  if (track.pt() > trackPtCut_) trackCounterEtaPt++;
379  }
380 
381  math::XYZPoint v1(vx,vy, vz);
382  double dz= track.dz(v1);
383  double dzsigma = sqrt(track.dzError()*track.dzError()+vzError*vzError);
384  double dxy= track.dxy(v1);
385  double dxysigma = sqrt(track.dxyError()*track.dxyError()+vxError*vyError);
386  if( track.quality(trackQuality_) &&
387  track.pt()>0.4 && fabs(track.eta())<2.4 &&
388  track.ptError()/track.pt()<0.1 && fabs(dz/dzsigma)<3.0 &&
389  fabs(dxy/dxysigma)<3.0){
390  nTracks++;
391  }
392  }
393 
394  creco->trackMultiplicity_ = nTracks;
395  creco->ntracksPtCut_ = trackCounter;
396  creco->ntracksEtaCut_ = trackCounterEta;
397  creco->ntracksEtaPtCut_ = trackCounterEtaPt;
398 
399  }else{
400  if(reuseAny_){
401  creco->trackMultiplicity_ = inputCentrality->Ntracks();
402  creco->ntracksPtCut_= inputCentrality->NtracksPtCut();
403  creco->ntracksEtaCut_= inputCentrality->NtracksEtaCut();
404  creco->ntracksEtaPtCut_= inputCentrality->NtracksEtaPtCut();
405  }
406  }
407 
408  if(producePixelTracks_){
409  Handle<TrackCollection> pixeltracks;
410  iEvent.getByToken(srcPixelTracks_,pixeltracks);
411  int nPixelTracks = pixeltracks->size();
412  creco->nPixelTracks_ = nPixelTracks;
413  }
414  else{
415  if(reuseAny_){
416  creco->nPixelTracks_ = inputCentrality->NpixelTracks();
417  }
418  }
419 
420  if(produceZDChits_){
421  creco->zdcSumPlus_ = 0;
422  creco->zdcSumMinus_ = 0;
423 
425  bool zdcAvailable = iEvent.getByToken(srcZDChits_,hits);
426  if(zdcAvailable){
427  for( size_t ihit = 0; ihit<hits->size(); ++ ihit){
428  const ZDCRecHit & rechit = (*hits)[ ihit ];
429  if(rechit.id().zside() > 0 ) {
430  if(lowGainZDC_){
431  creco->zdcSumPlus_ += rechit.lowGainEnergy();
432  }else{
433  creco->zdcSumPlus_ += rechit.energy();
434  }
435  }
436  if(rechit.id().zside() < 0) {
437  if(lowGainZDC_){
438  creco->zdcSumMinus_ += rechit.lowGainEnergy();
439  }else{
440  creco->zdcSumMinus_ += rechit.energy();
441  }
442  }
443  }
444  }else{
445  creco->zdcSumPlus_ = -9;
446  creco->zdcSumMinus_ = -9;
447  }
448  }else{
449  if(reuseAny_){
450  creco->zdcSumMinus_ = inputCentrality->zdcSumMinus();
451  creco->zdcSumPlus_ = inputCentrality->zdcSumPlus();
452  }
453  }
454 
455  iEvent.put(creco);
456 
457 }
458 
459 // ------------ method called once each job just before starting event loop ------------
460 void
462 {
463 }
464 
465 // ------------ method called once each job just after ending the event loop ------------
466 void
467 CentralityProducer::endJob()
468 {
469 }
470 
471 
472 //define this as a plug-in
474 
475 }
reco::TrackBase::TrackQuality trackQuality_
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
boost::transform_iterator< IterHelp, const_IdIter > const_iterator
const_iterator end(bool update=false) const
edm::EDGetTokenT< EcalRecHitCollection > srcEBhits_
virtual float pt() const
transverse momentum
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
edm::EDGetTokenT< Centrality > reuseTag_
TrackQuality
track quality
Definition: TrackBase.h:133
double dxyError() const
error on dxy
Definition: TrackBase.h:844
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:52
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
int zside() const
get the z-side of the cell (1/-1)
Definition: HcalZDCDetId.h:34
edm::EDGetTokenT< VertexCollection > srcVertex_
T eta() const
void beginJob()
Definition: Breakpoints.cc:15
edm::EDGetTokenT< CaloTowerCollection > srcTowers_
int iEvent
Definition: GenABIO.cc:230
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:699
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
float energy() const
Definition: CaloRecHit.h:17
virtual float eta() const
momentum pseudorapidity
T sqrt(T t)
Definition: SSEVec.h:48
double pt() const
track transverse momentum
Definition: TrackBase.h:669
double ptError() const
error on Pt (set to 1000 TeV if charge==0 for safety)
Definition: TrackBase.h:811
int ieta() const
get the cell ieta
Definition: HcalDetId.h:36
edm::EDGetTokenT< HFRecHitCollection > srcHFhits_
float energy() const
Definition: EcalRecHit.h:68
float lowGainEnergy() const
Definition: ZDCRecHit.h:21
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:657
double dzError() const
error on dz
Definition: TrackBase.h:862
Definition: DetId.h:18
bool isHF(int etabin, int depth)
DetId id() const
get the id
Definition: EcalRecHit.h:76
T const * product() const
Definition: Handle.h:81
tuple tracks
Definition: testEve_cfg.py:39
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::EDGetTokenT< ZDCRecHitCollection > srcZDChits_
ClusterRef cluster() const
Definition: SiPixelRecHit.h:49
const T & get() const
Definition: EventSetup.h:55
id_type detId() const
Definition: DetSetNew.h:83
edm::ESHandle< CaloGeometry > cGeo
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:551
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalPoint
point in global coordinate system
Definition: Point3D.h:17
const_iterator find(id_type i, bool update=false) const
edm::EDGetTokenT< SiPixelRecHitCollection > srcPixelhits_
int ietaAbs() const
Definition: CaloTower.h:170
iterator end()
Definition: DetSetNew.h:70
edm::EDGetTokenT< TrackCollection > srcPixelTracks_
edm::EDGetTokenT< TrackCollection > srcTracks_
HcalZDCDetId id() const
get the id
Definition: ZDCRecHit.h:19
edm::ESHandle< TrackerGeometry > tGeo
HcalDetId id() const
Definition: HFRecHit.h:23
DetId geographicalId() const
edm::EDGetTokenT< EcalRecHitCollection > srcEEhits_
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:639
virtual LocalPoint localPosition() const
const_iterator begin(bool update=false) const
Our base class.
Definition: SiPixelRecHit.h:23
iterator begin()
Definition: DetSetNew.h:67