CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EvtPlaneProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: EvtPlaneProducer
4 // Class: EvtPlaneProducer
5 //
13 //
14 // Original Author: Sergey Petrushanko
15 // Created: Fri Jul 11 10:05:00 2008
16 //
17 //
18 #define TRACKCOLLECTION 1
19 //#define RECOCHARGEDCANDIDATECOLLECTION 1
20 
21 
22 // system include files
23 #include <memory>
24 #include <iostream>
25 #include <time.h>
26 #include "TMath.h"
27 
28 // user include files
33 
37 
42 
46 
49 
56 
57 #include <cstdlib>
59 
60 using namespace std;
61 using namespace hi;
62 
63 //
64 // class decleration
65 //
66 
68 public:
69  explicit EvtPlaneProducer(const edm::ParameterSet&);
71 
72 private:
73  //edm::Service<TFileService> fs;
74  class GenPlane {
75  public:
76  GenPlane(string name,double etaminval1,double etamaxval1,double etaminval2,double etamaxval2,int orderval){
77  epname=name;
78  etamin1=etaminval1;
79  etamax1=etamaxval1;
80  etamin2=etaminval2;
81  etamax2=etamaxval2;
82  sumsin=0;
83  sumcos=0;
84  order = (double) orderval;
85  }
86  ~GenPlane(){;}
87  void addParticle(double w, double s, double c, double eta) {
88  if(w < 0.001) return;
89  if((eta>=etamin1 && eta<etamax1) ||
90  (etamin2!= etamax2 && eta>=etamin2 && eta<etamax2 )) {
91  sumsin+=w*s;
92  sumcos+=w*c;
93  //For tracking, w=1 and mult is the track multiplicity.
94  //For calorimetors, w is an energy that needs to be subsequently
95  //converted if an multiplicity if needed
96  mult+=w;
97  }
98  }
99 
100  double getAngle(double &ang, double &sv, double &cv){
101  ang = -10;
102  sv = 0;
103  cv = 0;
104  sv = sumsin;
105  cv = sumcos;
106  double q = sv*sv+cv*cv;
107  if(q>0) ang = atan2(sv,cv)/order;
108  return ang;
109  }
110  void reset() {
111  sumsin=0;
112  sumcos=0;
113  mult = 0;
114  }
115  private:
116  string epname;
117  double etamin1;
118  double etamax1;
119 
120  double etamin2;
121  double etamax2;
122  double sumsin;
123  double sumcos;
124  int mult;
125  double order;
126  };
127 
128 
130 
131  virtual void beginJob() override ;
132  virtual void produce(edm::Event&, const edm::EventSetup&) override;
133  virtual void endJob() override ;
134 
135  // ----------member data ---------------------------
139  bool useECAL_;
140  bool useHCAL_;
141  bool useTrack_;
143  double minet_;
144  double maxet_;
145  double minpt_;
146  double maxpt_;
147  double minvtx_;
148  double maxvtx_;
149  double dzerr_;
150  double chi2_;
151 
153 };
154 
155 //
156 // constants, enums and typedefs
157 //
158 
159 
160 //
161 // static data member definitions
162 //
163 
164 //
165 // constructors and destructor
166 //
168 {
169 
170 
171  //register your products
172  vtxCollection_ = iConfig.getParameter<edm::InputTag>("vtxCollection_");
173  caloCollection_ = iConfig.getParameter<edm::InputTag>("caloCollection_");
174  trackCollection_ = iConfig.getParameter<edm::InputTag>("trackCollection_");
175  useECAL_ = iConfig.getUntrackedParameter<bool>("useECAL_",true);
176  useHCAL_ = iConfig.getUntrackedParameter<bool>("useHCAL_",true);
177  useTrack_ = iConfig.getUntrackedParameter<bool>("useTrack",true);
178  useTrackPtWeight_ = iConfig.getUntrackedParameter<bool>("useTrackPtWeight_",true);
179  minet_ = iConfig.getUntrackedParameter<double>("minet_",0.2);
180  maxet_ = iConfig.getUntrackedParameter<double>("maxet_",500.);
181  minpt_ = iConfig.getUntrackedParameter<double>("minpt_",0.3);
182  maxpt_ = iConfig.getUntrackedParameter<double>("maxpt_",2.5);
183  minvtx_ = iConfig.getUntrackedParameter<double>("minvtx_",-50.);
184  maxvtx_ = iConfig.getUntrackedParameter<double>("maxvtx_",50.);
185  dzerr_ = iConfig.getUntrackedParameter<double>("dzerr_",10.);
186  chi2_ = iConfig.getUntrackedParameter<double>("chi2_",40.);
187 
188  storeNames_ = 1;
189  produces<reco::EvtPlaneCollection>("recoLevel");
190  for(int i = 0; i<NumEPNames; i++ ) {
191  rp[i] = new GenPlane(EPNames[i].data(),EPEtaMin1[i],EPEtaMax1[i],EPEtaMin2[i],EPEtaMax2[i],EPOrder[i]);
192  }
193 }
194 
195 
197 {
198 
199  // do anything here that needs to be done at desctruction time
200  // (e.g. close files, deallocate resources etc.)
201 
202 }
203 
204 
205 //
206 // member functions
207 //
208 
209 // ------------ method called to produce the data ------------
210 void
212 {
213  using namespace edm;
214  using namespace reco;
215 
216  int vs_sell;
217  float vzr_sell;
218  //
219  //Get Vertex
220  //
221  edm::Handle<reco::VertexCollection> vertexCollection3;
222  iEvent.getByLabel(vtxCollection_,vertexCollection3);
223  const reco::VertexCollection * vertices3 = vertexCollection3.product();
224  vs_sell = vertices3->size();
225  if(vs_sell>0) {
226  vzr_sell = vertices3->begin()->z();
227  } else
228  vzr_sell = -999.9;
229  //
230  for(int i = 0; i<NumEPNames; i++) rp[i]->reset();
231  if(vzr_sell>minvtx_ && vzr_sell<maxvtx_) {
232  //calorimetry part
233 
234  double tower_eta, tower_phi;
235  double tower_energyet, tower_energyet_e, tower_energyet_h;
236  double s1, s2, s13,s23,s14,s24,s15,s25,s16,s26;
237  Handle<CaloTowerCollection> calotower;
238  iEvent.getByLabel(caloCollection_,calotower);
239 
240  if(calotower.isValid()){
241 
242  for (CaloTowerCollection::const_iterator j = calotower->begin();j !=calotower->end(); j++) {
243  tower_eta = j->eta();
244  tower_phi = j->phi();
245  tower_energyet_e = j->emEt();
246  tower_energyet_h = j->hadEt();
247  tower_energyet = tower_energyet_e + tower_energyet_h;
248 
249  s1 = sin(2.*tower_phi);
250  s2 = cos(2.*tower_phi);
251  s13 = sin(3.*tower_phi);
252  s23 = cos(3.*tower_phi);
253  s14 = sin(4.*tower_phi);
254  s24 = cos(4.*tower_phi);
255  s15 = sin(5.*tower_phi);
256  s25 = cos(5.*tower_phi);
257  s16 = sin(6.*tower_phi);
258  s26 = cos(6.*tower_phi);
259 
260  if(tower_energyet<minet_) continue;
261  if(tower_energyet>maxet_) continue;;
262 
263  rp[etEcal ]->addParticle (tower_energyet_e, s1, s2, tower_eta);
264  rp[etEcalP ]->addParticle (tower_energyet_e, s1, s2, tower_eta);
265  rp[etEcalM ]->addParticle (tower_energyet_e, s1, s2, tower_eta);
266  rp[etHcal ]->addParticle (tower_energyet_h, s1, s2, tower_eta);
267  rp[etHcalP ]->addParticle (tower_energyet_h, s1, s2, tower_eta);
268  rp[etHcalM ]->addParticle (tower_energyet_h, s1, s2, tower_eta);
269  rp[etHF ]->addParticle (tower_energyet, s1, s2, tower_eta);
270  rp[etHFp ]->addParticle (tower_energyet, s1, s2, tower_eta);
271  rp[etHFm ]->addParticle (tower_energyet, s1, s2, tower_eta);
272  rp[etCaloHFP ]->addParticle (tower_energyet, s1, s2, tower_eta);
273  rp[etCaloHFM ]->addParticle (tower_energyet, s1, s2, tower_eta);
274 
275  rp[etHF3 ]->addParticle (tower_energyet, s13, s23, tower_eta);
276  rp[etHFp3 ]->addParticle (tower_energyet, s13, s23, tower_eta);
277  rp[etHFm3 ]->addParticle (tower_energyet, s13, s23, tower_eta);
278 
279  rp[etHF4 ]->addParticle (tower_energyet, s14, s24, tower_eta);
280  rp[etHFp4 ]->addParticle (tower_energyet, s14, s24, tower_eta);
281  rp[etHFm4 ]->addParticle (tower_energyet, s14, s24, tower_eta);
282 
283  rp[etHF5 ]->addParticle (tower_energyet, s15, s25, tower_eta);
284  rp[etHFp5 ]->addParticle (tower_energyet, s15, s25, tower_eta);
285  rp[etHFm5 ]->addParticle (tower_energyet, s15, s25, tower_eta);
286 
287  rp[etHF6 ]->addParticle (tower_energyet, s16, s26, tower_eta);
288  rp[etHFp6 ]->addParticle (tower_energyet, s16, s26, tower_eta);
289  rp[etHFm6 ]->addParticle (tower_energyet, s16, s26, tower_eta);
290  }
291  }
292  //Tracking part
293 
294  double track_eta;
295  double track_phi;
296  double track_pt;
297 #ifdef TRACKCOLLECTION
299  iEvent.getByLabel(trackCollection_, tracks);
300 
301  if(tracks.isValid()){
302  for(reco::TrackCollection::const_iterator j = tracks->begin(); j != tracks->end(); j++){
303 
304  //Find possible collections with command line: edmDumpEventContent *.root
305 #endif
306 #ifdef RECOCHARGEDCANDIDATECOLLECTION
308  iEvent.getByLabel("allMergedPtSplit12Tracks",trackCollection);
309  // iEvent.getByLabel("hiGoodMergedTracks",trackCollection);
310  if(trackCollection.isValid()){
311 
312  const reco::RecoChargedCandidateCollection * tracks = trackCollection.product();
313  for(reco::RecoChargedCandidateCollection::const_iterator j = tracks->begin(); j != tracks->end(); j++){
314 #endif
316  iEvent.getByLabel(vtxCollection_, vertex);
317 
318 // find the vertex point and error
319 
320  math::XYZPoint vtxPoint(0.0,0.0,0.0);
321  double vzErr =0.0, vxErr=0.0, vyErr=0.0;
322  if(vertex->size()>0) {
323  vtxPoint=vertex->begin()->position();
324  vzErr=vertex->begin()->zError();
325  vxErr=vertex->begin()->xError();
326  vyErr=vertex->begin()->yError();
327  }
328  bool accepted = true;
329  bool isPixel = false;
330  // determine if the track is a pixel track
331  if ( j->numberOfValidHits() < 7 ) isPixel = true;
332 
333  // determine the vertex significance
334  double d0=0.0, dz=0.0, d0sigma=0.0, dzsigma=0.0;
335  d0 = -1.*j->dxy(vtxPoint);
336  dz = j->dz(vtxPoint);
337  d0sigma = sqrt(j->d0Error()*j->d0Error()+vxErr*vyErr);
338  dzsigma = sqrt(j->dzError()*j->dzError()+vzErr*vzErr);
339 
340  // cuts for pixel tracks
341  if( isPixel )
342  {
343  // dz significance cut
344  if ( fabs(dz/dzsigma) > dzerr_ ) accepted = false;
345 
346  // chi2/ndof cut
347  if ( j->normalizedChi2() > chi2_ ) accepted = false;
348  }
349 
350  // cuts for full tracks
351  if ( ! isPixel)
352  {
353  // dz and d0 significance cuts
354  if ( fabs(dz/dzsigma) > 3 ) accepted = false;
355  if ( fabs(d0/d0sigma) > 3 ) accepted = false;
356 
357  // pt resolution cut
358  if ( j->ptError()/j->pt() > 0.05 ) accepted = false;
359 
360  // number of valid hits cut
361  if ( j->numberOfValidHits() < 12 ) accepted = false;
362  }
363  if( accepted ) {
364  track_eta = j->eta();
365  track_phi = j->phi();
366  track_pt = j->pt();
367  double s =sin(2*track_phi);
368  double c =cos(2*track_phi);
369  double s3 =sin(3*track_phi);
370  double c3 =cos(3*track_phi);
371  double s4 =sin(4*track_phi);
372  double c4 =cos(4*track_phi);
373  double s5 =sin(5*track_phi);
374  double c5 =cos(5*track_phi);
375  double s6 =sin(6*track_phi);
376  double c6 =cos(6*track_phi);
377  double w = 1;
378  if(useTrackPtWeight_) {
379  w = track_pt;
380  if(w>2.5) w=2.0; //v2 starts decreasing above ~2.5 GeV/c
381  }
382  if(track_pt<minpt_) continue;
383  if(track_pt>maxpt_) continue;
384  rp[EvtPlaneFromTracksMidEta]->addParticle(w,s,c,track_eta);
385  rp[EvtPTracksPosEtaGap]->addParticle(w,s,c,track_eta);
386  rp[EvtPTracksNegEtaGap]->addParticle(w,s,c,track_eta);
387  rp[EPTracksMid3]->addParticle(w,s3,c3,track_eta);
388  rp[EPTracksPos3]->addParticle(w,s3,c3,track_eta);
389  rp[EPTracksNeg3]->addParticle(w,s3,c3,track_eta);
390  rp[EPTracksMid4]->addParticle(w,s4,c4,track_eta);
391  rp[EPTracksPos4]->addParticle(w,s4,c4,track_eta);
392  rp[EPTracksNeg4]->addParticle(w,s4,c4,track_eta);
393  rp[EPTracksMid5]->addParticle(w,s5,c5,track_eta);
394  rp[EPTracksPos5]->addParticle(w,s5,c5,track_eta);
395  rp[EPTracksNeg5]->addParticle(w,s5,c5,track_eta);
396  rp[EPTracksMid6]->addParticle(w,s6,c6,track_eta);
397  rp[EPTracksPos6]->addParticle(w,s6,c6,track_eta);
398  rp[EPTracksNeg6]->addParticle(w,s6,c6,track_eta);
399 
400  }
401 
402  }
403  }
404  }
405  std::auto_ptr<EvtPlaneCollection> evtplaneOutput(new EvtPlaneCollection);
406 
407  EvtPlane *ep[NumEPNames];
408 
409  double ang=-10;
410  double sv = 0;
411  double cv = 0;
412 
413  for(int i = 0; i<NumEPNames; i++) {
414  rp[i]->getAngle(ang,sv,cv);
415  if(storeNames_) ep[i] = new EvtPlane(ang,sv,cv,EPNames[i]);
416  else ep[i] = new EvtPlane(ang,sv,cv,"");
417  }
418  if(useTrack_) {
419  for(int i = 0; i<9; i++) {
420  evtplaneOutput->push_back(*ep[i]);
421  }
422  }
423  for(int i = 9; i<NumEPNames; i++) {
424  if(useECAL_ && !useHCAL_) {
425  if(EPNames[i].rfind("Ecal")!=string::npos) {
426  evtplaneOutput->push_back(*ep[i]);
427  }
428  } else if (useHCAL_ && !useECAL_) {
429  if(EPNames[i].rfind("Hcal")!=string::npos) {
430  evtplaneOutput->push_back(*ep[i]);
431  }
432  }else if (useECAL_ && useHCAL_) {
433  evtplaneOutput->push_back(*ep[i]);
434  }
435  }
436 
437  iEvent.put(evtplaneOutput, "recoLevel");
438  // storeNames_ = 0;
439 }
440 
441  // ------------ method called once each job just before starting event loop ------------
442 void
444 {
445 }
446 
447 // ------------ method called once each job just after ending the event loop ------------
448 void
450 }
451 
452 //define this as a plug-in
double getAngle(double &ang, double &sv, double &cv)
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const int EPOrder[]
int i
Definition: DBlmapReader.cc:9
const double EPEtaMin2[]
GenPlane(string name, double etaminval1, double etamaxval1, double etaminval2, double etamaxval2, int orderval)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void beginJob() override
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< CaloTower >::const_iterator const_iterator
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
T eta() const
const std::string EPNames[]
virtual void produce(edm::Event &, const edm::EventSetup &) override
void beginJob()
Definition: Breakpoints.cc:15
tuple s2
Definition: indexGen.py:106
const double EPEtaMin1[]
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
T sqrt(T t)
Definition: SSEVec.h:48
virtual void endJob() override
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::vector< EvtPlane > EvtPlaneCollection
Definition: EvtPlane.h:33
int j
Definition: DBlmapReader.cc:9
EvtPlaneProducer(const edm::ParameterSet &)
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
dictionary cv
Definition: cuy.py:362
std::vector< RecoChargedCandidate > RecoChargedCandidateCollection
collectin of RecoChargedCandidate objects
void addParticle(double w, double s, double c, double eta)
tuple tracks
Definition: testEve_cfg.py:39
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
T const * product() const
Definition: Handle.h:81
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
edm::InputTag caloCollection_
edm::InputTag vtxCollection_
const double EPEtaMax2[]
static const int NumEPNames
T w() const
edm::InputTag trackCollection_
void reset(double vett[256])
Definition: TPedValues.cc:11
const double EPEtaMax1[]