CMS 3D CMS Logo

DAClusterizerInZT_vect.h
Go to the documentation of this file.
1 #ifndef DAClusterizerInZT_vect_h
2 #define DAClusterizerInZT_vect_h
3 
15 #include <vector>
19 
20 #include <memory>
21 
23 
24 public:
25 
26  // Internal data structure to
27  struct track_t {
28 
29  void addItem( double new_z, double new_t, double new_dz2, double new_dt2, const reco::TransientTrack* new_tt, double new_pi )
30  {
31  z.push_back( new_z );
32  t.push_back( new_t );
33  dz2.push_back( new_dz2 );
34  dt2.push_back( new_dt2 );
35  tt.push_back( new_tt );
36 
37  pi.push_back( new_pi ); // track weight
38  Z_sum.push_back( 1.0 ); // Z[i] for DA clustering, initial value as done in ::fill
39  }
40 
41  unsigned int getSize() const
42  {
43  return z.size();
44  }
45 
46  // has to be called everytime the items are modified
47  void extractRaw()
48  {
49  z_ = &z.front();
50  t_ = &t.front();
51  dz2_ = &dz2.front();
52  dt2_ = &dt2.front();
53  Z_sum_ = &Z_sum.front();
54  pi_ = &pi.front();
55  }
56 
57  double * z_; // z-coordinate at point of closest approach to the beamline
58  double * t_; // t-coordinate at point of closest approach to the beamline
59  double * pi_; // track weight
60 
61  double * dz2_; // square of the error of z(pca)
62  double * dt2_; // square of the error of t(pca)
63  double * Z_sum_; // Z[i] for DA clustering
64 
65  std::vector<double> z; // z-coordinate at point of closest approach to the beamline
66  std::vector<double> t; // t-coordinate at point of closest approach to the beamline
67  std::vector<double> dz2; // square of the error of z(pca)
68  std::vector<double> dt2; // square of the error of t(pca)
69  std::vector<double> Z_sum; // Z[i] for DA clustering
70  std::vector<double> pi; // track weight
71  std::vector< const reco::TransientTrack* > tt; // a pointer to the Transient Track
72  };
73 
74  struct vertex_t {
75 
76  void addItem( double new_z, double new_t, double new_pk )
77  {
78  z.push_back( new_z);
79  t.push_back( new_t);
80  pk.push_back( new_pk);
81 
82  ei_cache.push_back( 0.0 );
83  ei.push_back( 0.0 );
84  swz.push_back( 0.0);
85  swt.push_back( 0.0);
86  se.push_back( 0.0);
87  nuz.push_back(0.0);
88  nut.push_back(0.0);
89  szz.push_back(0.0);
90  stt.push_back(0.0);
91  szt.push_back(0.0);
92 
93  extractRaw();
94  }
95 
96  unsigned int getSize() const
97  {
98  return z.size();
99  }
100 
101  // has to be called everytime the items are modified
102  void extractRaw()
103  {
104  z_ = &z.front();
105  t_ = &t.front();
106  pk_ = &pk.front();
107 
108  ei_ = &ei.front();
109  swz_ = &swz.front();
110  swt_ = &swt.front();
111  se_ = &se.front();
112  nuz_ = &nuz.front();
113  nut_ = &nut.front();
114  szz_ = &szz.front();
115  stt_ = &stt.front();
116  szt_ = &szt.front();
117 
118  ei_cache_ = &ei_cache.front();
119 
120  }
121 
122  void insertItem( unsigned int i, double new_z, double new_t, double new_pk )
123  {
124  z.insert(z.begin() + i, new_z);
125  t.insert(t.begin() + i, new_t);
126  pk.insert(pk.begin() + i, new_pk);
127 
128  ei_cache.insert(ei_cache.begin() + i, 0.0 );
129  ei.insert( ei.begin() + i, 0.0 );
130  swz.insert(swz.begin() + i, 0.0 );
131  swt.insert(swt.begin() + i, 0.0 );
132  se.insert( se.begin() + i, 0.0 );
133 
134  nuz.insert(nuz.begin() +i, 0.0 );
135  nut.insert(nut.begin() +i, 0.0 );
136  szz.insert(szz.begin() + i, 0.0 );
137  stt.insert(stt.begin() + i, 0.0 );
138  szt.insert(szt.begin() + i, 0.0 );
139  extractRaw();
140  }
141 
142  void removeItem( unsigned int i )
143  {
144  z.erase( z.begin() + i );
145  t.erase( t.begin() + i );
146  pk.erase( pk.begin() + i );
147 
148  ei_cache.erase( ei_cache.begin() + i);
149  ei.erase( ei.begin() + i);
150  swz.erase( swz.begin() + i);
151  swt.erase( swt.begin() + i);
152  se.erase(se.begin() + i);
153 
154  nuz.erase(nuz.begin() + i);
155  nut.erase(nut.begin() + i);
156  szz.erase(szz.begin() + i);
157  stt.erase(stt.begin() + i);
158  szt.erase(szt.begin() + i);
159 
160  extractRaw();
161  }
162 
163 
164  unsigned int insertOrdered( double z, double t, double pk){
165  // insert a new cluster according to it's z-position, return the index at which it was inserted
166 
167  unsigned int k = 0;
168  for( ; k < getSize(); k++){
169  if (z < z_[k]) break;
170  }
171  insertItem(k ,z, t, pk);
172  return k;
173  }
174 
175 
176  void debugOut()
177  {
178  std::cout << "vertex_t size: " << getSize() << std::endl;
179 
180  for ( unsigned int i =0; i < getSize(); ++ i)
181  {
182  std::cout << " z = " << z_[i] << " t = " << t_[i] << " pk = " << pk_[i] << std::endl;
183  }
184  }
185 
186  std::vector<double> z; // z coordinate
187  std::vector<double> t; // t coordinate
188  std::vector<double> pk; // vertex weight for "constrained" clustering
189 
190  double * z_;
191  double * t_;
192  double * pk_;
193 
194  double * ei_cache_;
195  double * ei_;
196  double * swz_;
197  double * swt_;
198  double * se_;
199  double * szz_;
200  double * stt_;
201  double * szt_;
202  double * nuz_;
203  double * nut_;
204 
205  // --- temporary numbers, used during update
206  std::vector<double> ei_cache;
207  std::vector<double> ei;
208  std::vector<double> swz;
209  std::vector<double> swt;
210  std::vector<double> se;
211  std::vector<double> nuz;
212  std::vector<double> nut;
213  std::vector<double> szz;
214  std::vector<double> stt;
215  std::vector<double> szt;
216  };
217 
219 
220  std::vector<std::vector<reco::TransientTrack> >
221  clusterize(const std::vector<reco::TransientTrack> & tracks) const override;
222 
223  std::vector<TransientVertex>
224  vertices(const std::vector<reco::TransientTrack> & tracks,
225  const int verbosity = 0) const ;
226 
227  track_t fill(const std::vector<reco::TransientTrack> & tracks) const;
228 
229  double update(double beta, track_t & gtracks,
230  vertex_t & gvertices, bool useRho0, const double & rho0) const;
231 
232  void dump(const double beta, const vertex_t & y,
233  const track_t & tks, const int verbosity = 0) const;
234  void zorder(vertex_t & y)const;
235  bool find_nearest(double z, double t, vertex_t & y, unsigned int & k_min, double dz, double dt)const;
236  bool merge(vertex_t & y, double & beta)const;
237  bool purge(vertex_t &, track_t &, double &,
238  const double) const;
239  void splitAll( vertex_t & y) const;
240  bool split(const double beta, track_t &t, vertex_t & y, double threshold = 1. ) const;
241 
242  double beta0(const double betamax, track_t const & tks, vertex_t const & y) const;
243 
244  double get_Tc(const vertex_t & y, int k) const;
245 
246 private:
247  bool verbose_;
248  double zdumpcenter_;
249  double zdumpwidth_;
250 
251  double vertexSize_;
255  double betamax_;
256  double betastop_;
257  double dzCutOff_;
258  double d0CutOff_;
259  double dtCutOff_;
260  bool useTc_;
261 
264  double zmerge_;
265  double tmerge_;
266  double betapurge_;
267 
268 };
269 
270 
271 //#ifndef DAClusterizerInZT_new_h
272 #endif
const double beta
float dt
Definition: AMPTWrapper.h:126
unsigned int insertOrdered(double z, double t, double pk)
std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const override
bool split(const double beta, track_t &t, vertex_t &y, double threshold=1.) const
void splitAll(vertex_t &y) const
DAClusterizerInZT_vect(const edm::ParameterSet &conf)
bool merge(vertex_t &y, double &beta) const
void addItem(double new_z, double new_t, double new_dz2, double new_dt2, const reco::TransientTrack *new_tt, double new_pi)
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks, const int verbosity=0) const
void zorder(vertex_t &y) const
void insertItem(unsigned int i, double new_z, double new_t, double new_pk)
double get_Tc(const vertex_t &y, int k) const
bool purge(vertex_t &, track_t &, double &, const double) const
bool find_nearest(double z, double t, vertex_t &y, unsigned int &k_min, double dz, double dt) const
double update(double beta, track_t &gtracks, vertex_t &gvertices, bool useRho0, const double &rho0) const
int k[5][pyjets_maxn]
void addItem(double new_z, double new_t, double new_pk)
track_t fill(const std::vector< reco::TransientTrack > &tracks) const
double beta0(const double betamax, track_t const &tks, vertex_t const &y) const
std::vector< const reco::TransientTrack * > tt
void dump(const double beta, const vertex_t &y, const track_t &tks, const int verbosity=0) const