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