CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTTrack.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <sstream>
3 #include <istream>
4 #include <fstream>
5 #include <iomanip>
6 #include <string>
7 #include <cmath>
8 #include <functional>
9 #include <stdlib.h>
10 #include <string.h>
11 
14 
16  evtCounter=0;
17 
18  //set parameter defaults
19  _Monte=false;
20  _Debug=false;
21 }
22 
23 /* Setup the analysis to put the branch-variables into the tree. */
24 void HLTTrack::setup(const edm::ParameterSet& pSet, TTree* HltTree) {
25 
26  edm::ParameterSet myEmParams = pSet.getParameter<edm::ParameterSet>("RunParameters") ;
27  std::vector<std::string> parameterNames = myEmParams.getParameterNames() ;
28 
29  for ( std::vector<std::string>::iterator iParam = parameterNames.begin();
30  iParam != parameterNames.end(); iParam++ ){
31  if ( (*iParam) == "Monte" ) _Monte = myEmParams.getParameter<bool>( *iParam );
32  else if ( (*iParam) == "Debug" ) _Debug = myEmParams.getParameter<bool>( *iParam );
33  }
34 
35  //common
36  const int kMaxTrackL3 = 10000;
37  //isoPixel
38  isopixeltrackL3pt = new float[kMaxTrackL3];
39  isopixeltrackL3eta = new float[kMaxTrackL3];
40  isopixeltrackL3phi = new float[kMaxTrackL3];
41  isopixeltrackL3maxptpxl = new float[kMaxTrackL3];
42  isopixeltrackL3energy = new float[kMaxTrackL3];
43  isopixeltrackL2pt = new float[kMaxTrackL3];
44  isopixeltrackL2eta = new float[kMaxTrackL3];
45  isopixeltrackL2dXY = new float[kMaxTrackL3];
46 
47  //minBiasPixel
48  pixeltracksL3pt = new float[kMaxTrackL3];
49  pixeltracksL3eta = new float[kMaxTrackL3];
50  pixeltracksL3phi = new float[kMaxTrackL3];
51  pixeltracksL3vz = new float[kMaxTrackL3];
52 
53  // Track-specific branches of the tree
54  //isoPixel
55  HltTree->Branch("NohIsoPixelTrackL3",&nisopixeltrackL3,"NohIsoPixelTrackL3/I");
56  HltTree->Branch("ohIsoPixelTrackL3Pt",isopixeltrackL3pt,"ohIsoPixelTrackL3Pt[NohIsoPixelTrackL3]/F");
57  HltTree->Branch("ohIsoPixelTrackL3Eta",isopixeltrackL3eta,"ohIsoPixelTrackL3Eta[NohIsoPixelTrackL3]/F");
58  HltTree->Branch("ohIsoPixelTrackL3Phi",isopixeltrackL3phi,"ohIsoPixelTrackL3Phi[NohIsoPixelTrackL3]/F");
59  HltTree->Branch("ohIsoPixelTrackL3MaxPtPxl",isopixeltrackL3maxptpxl,"ohIsoPixelTrackL3MaxPtPxl[NohIsoPixelTrackL3]/F");
60  HltTree->Branch("ohIsoPixelTrackL3Energy",isopixeltrackL3energy,"ohIsoPixelTrackL3Energy[NohIsoPixelTrackL3]/F");
61  HltTree->Branch("ohIsoPixelTrackL2pt",isopixeltrackL2pt,"ohIsoPixelTrackL2pt[NohIsoPixelTrackL3]/F");
62  HltTree->Branch("ohIsoPixelTrackL2eta",isopixeltrackL2eta,"ohIsoPixelTrackL2eta[NohIsoPixelTrackL3]/F");
63  HltTree->Branch("ohIsoPixelTrackL2dXY",isopixeltrackL2dXY,"ohIsoPixelTrackL2dXY[NohIsoPixelTrackL3]/F");
64 
65  //minBiasPixel
66  HltTree->Branch("NohPixelTracksL3",&npixeltracksL3,"NohPixelTracksL3/I");
67  HltTree->Branch("ohPixelTracksL3Pt",pixeltracksL3pt,"ohPixelTracksL3Pt[NohPixelTracksL3]/F");
68  HltTree->Branch("ohPixelTracksL3Eta",pixeltracksL3eta,"ohPixelTracksL3Eta[NohPixelTracksL3]/F");
69  HltTree->Branch("ohPixelTracksL3Phi",pixeltracksL3phi,"ohPixelTracksL3Phi[NohPixelTracksL3]/F");
70  HltTree->Branch("ohPixelTracksL3Vz",pixeltracksL3vz,"ohPixelTracksL3Vz[NohPixelTracksL3]/F");
71 }
72 
73 /* **Analyze the event** */
79  TTree* HltTree) {
80 
81  //isoPixel
82  if (IsoPixelTrackL3.isValid()) {
83  // Ref to Candidate object to be recorded in filter object
85 
86  nisopixeltrackL3 = IsoPixelTrackL3->size();
87 
88  for (unsigned int i=0; i<IsoPixelTrackL3->size(); i++)
89  {
90  candref = edm::Ref<reco::IsolatedPixelTrackCandidateCollection>(IsoPixelTrackL3, i);
91 
92  isopixeltrackL3maxptpxl[i] = candref->maxPtPxl();
93  isopixeltrackL3pt[i] = candref->pt();
94  isopixeltrackL3eta[i] = candref->track()->eta();
95  isopixeltrackL3phi[i] = candref->phi();
96  isopixeltrackL3energy[i] = (candref->pt())*cosh(candref->track()->eta());
97  }
98  }
99  else {nisopixeltrackL3 = 0;}
100 
101  //minBiasPixel
102  if (PixelTracksL3.isValid()) {
103  // Ref to Candidate object to be recorded in filter object
105 
106  npixeltracksL3 = PixelTracksL3->size();
107 
108  for (unsigned int i=0; i<PixelTracksL3->size(); i++)
109  {
110  candref = edm::Ref<reco::RecoChargedCandidateCollection>(PixelTracksL3, i);
111 
112  pixeltracksL3pt[i] = candref->pt();
113  pixeltracksL3eta[i] = candref->eta();
114  pixeltracksL3phi[i] = candref->phi();
115  pixeltracksL3vz[i] = candref->vz();
116 
117  if (IsoPixelTrackL2.isValid()) {
118  double minDR=100;
121  for (unsigned int j=0; j<IsoPixelTrackL2->size(); j++)
122  {
123  candrefl2 = edm::Ref<reco::IsolatedPixelTrackCandidateCollection>(IsoPixelTrackL2, j);
124  double drL3L2 = deltaR(candrefl2->eta(), candrefl2->phi(),candref->eta(), candref->phi());
125  if (drL3L2<minDR)
126  {
127  candrefl2matched=candrefl2;
128  minDR=drL3L2;
129  }
130  }
131  if (candrefl2matched.isNonnull())
132  {
133  isopixeltrackL2pt[i]=candrefl2matched->pt();
134  isopixeltrackL2eta[i]=candrefl2matched->eta();
135  if (pixelVertices.isValid())
136  {
137  double minDZ=100;
139  edm::Ref<reco::VertexCollection> vertrefMatched;
140  for (unsigned int k=0; k<pixelVertices->size(); k++)
141  {
143  double dz=fabs(candrefl2matched->track()->dz(vertref->position()));
144  if (dz<minDZ)
145  {
146  minDZ=dz;
147  vertrefMatched=vertref;
148  }
149  }
150  if (vertrefMatched.isNonnull()) isopixeltrackL2dXY[i] = candrefl2matched->track()->dxy(vertref->position());
151  else isopixeltrackL2dXY[i]=0;
152  }
153  }
154  }
155  }
156  }
157  else {npixeltracksL3 = 0;}
158 
159 
161 
162 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
float * isopixeltrackL3pt
Definition: HLTTrack.h:52
bool _Monte
Definition: HLTTrack.h:59
void analyze(const edm::Handle< reco::IsolatedPixelTrackCandidateCollection > &IsoPixelTrackL3, const edm::Handle< reco::IsolatedPixelTrackCandidateCollection > &IsoPixelTrackL2, const edm::Handle< reco::VertexCollection > &pixelVertices, const edm::Handle< reco::RecoChargedCandidateCollection > &PixelTracksL3, TTree *tree)
Definition: HLTTrack.cc:74
float * pixeltracksL3eta
Definition: HLTTrack.h:55
float * pixeltracksL3pt
Definition: HLTTrack.h:55
int nisopixeltrackL3
Definition: HLTTrack.h:53
float * isopixeltrackL3maxptpxl
Definition: HLTTrack.h:52
float * isopixeltrackL2pt
Definition: HLTTrack.h:52
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:249
float * isopixeltrackL2eta
Definition: HLTTrack.h:52
float * isopixeltrackL3eta
Definition: HLTTrack.h:52
float * isopixeltrackL2dXY
Definition: HLTTrack.h:52
int j
Definition: DBlmapReader.cc:9
void setup(const edm::ParameterSet &pSet, TTree *tree)
Definition: HLTTrack.cc:24
std::vector< std::string > getParameterNames() const
bool isValid() const
Definition: HandleBase.h:76
int k[5][pyjets_maxn]
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
int npixeltracksL3
Definition: HLTTrack.h:56
float * pixeltracksL3vz
Definition: HLTTrack.h:55
float * isopixeltrackL3phi
Definition: HLTTrack.h:52
float * pixeltracksL3phi
Definition: HLTTrack.h:55
float * isopixeltrackL3energy
Definition: HLTTrack.h:52
HLTTrack()
Definition: HLTTrack.cc:15
int evtCounter
Definition: HLTTrack.h:61
bool _Debug
Definition: HLTTrack.h:59