CMS 3D CMS Logo

DAClusterizerInZT_vect.h
Go to the documentation of this file.
1 #ifndef RecoVertex_PrimaryVertexProducer_DAClusterizerInZT_vect_h
2 #define RecoVertex_PrimaryVertexProducer_DAClusterizerInZT_vect_h
3 
17 #include <vector>
21 
22 #include <memory>
23 
24 //#define USEVTXDT2
25 
27 public:
28  // internal data structure for tracks
29  struct track_t {
30  std::vector<double> zpca_vec; // z-coordinate at point of closest approach to the beamline
31  std::vector<double> tpca_vec; // t-coordinate at point of closest approach to the beamline
32  std::vector<double> dz2_vec; // square of the error of z(pca)
33  std::vector<double> dt2_vec; // square of the error of t(pca)
34  std::vector<double> sum_Z_vec; // track contribution to the partition function, Z
35  std::vector<double> tkwt_vec; // track weight, close to 1.0 for most tracks
36  std::vector<unsigned int> kmin; // index of the first cluster within zrange
37  std::vector<unsigned int> kmax; // 1 + index of the last cluster within zrange
38  std::vector<const reco::TransientTrack *> tt; // a pointer to the Transient Track
39 
40  double osumtkwt; // 1. / (sum of all track weights)
41 
42  void addItem(double new_zpca,
43  double new_tpca,
44  double new_dz2,
45  double new_dt2,
46  const reco::TransientTrack *new_tt,
47  double new_tkwt) {
48  zpca_vec.push_back(new_zpca);
49  tpca_vec.push_back(new_tpca);
50  dz2_vec.push_back(new_dz2);
51  dt2_vec.push_back(new_dt2);
52  tt.push_back(new_tt);
53  tkwt_vec.push_back(new_tkwt);
54  sum_Z_vec.push_back(1.0);
55  kmin.push_back(0);
56  kmax.push_back(0);
57  }
58 
59  void insertItem(unsigned int i,
60  double new_zpca,
61  double new_tpca,
62  double new_dz2,
63  double new_dt2,
64  const reco::TransientTrack *new_tt,
65  double new_tkwt) {
66  zpca_vec.insert(zpca_vec.begin() + i, new_zpca);
67  tpca_vec.insert(tpca_vec.begin() + i, new_tpca);
68  dz2_vec.insert(dz2_vec.begin() + i, new_dz2);
69  dt2_vec.insert(dt2_vec.begin() + i, new_dt2);
70  tt.insert(tt.begin() + i, new_tt);
71  tkwt_vec.insert(tkwt_vec.begin() + i, new_tkwt);
72  sum_Z_vec.insert(sum_Z_vec.begin() + i, 1.0);
73  kmin.insert(kmin.begin() + i, 0);
74  kmax.insert(kmax.begin() + i, 0);
75  }
76 
77  unsigned int getSize() const { return zpca_vec.size(); }
78 
79  // has to be called everytime the items are modified
80  void extractRaw() {
81  zpca = &zpca_vec.front();
82  tpca = &tpca_vec.front();
83  dz2 = &dz2_vec.front();
84  dt2 = &dt2_vec.front();
85  tkwt = &tkwt_vec.front();
86  sum_Z = &sum_Z_vec.front();
87  }
88 
89  // pointers to the first element of vectors, needed for vectorized code
90  double *__restrict__ zpca;
91  double *__restrict__ tpca;
92  double *__restrict__ dz2;
93  double *__restrict__ dt2;
94  double *__restrict__ tkwt;
95  double *__restrict__ sum_Z;
96  };
97 
98  // internal data structure for clusters
99  struct vertex_t {
100  std::vector<double> zvtx_vec; // z coordinate
101  std::vector<double> tvtx_vec; // t coordinate
102  std::vector<double> rho_vec; // vertex "mass" for mass-constrained clustering
103 #ifdef USEVTXDT2
104  std::vector<double> dt2_vec; // only used with vertex time uncertainties
105  std::vector<double> sumw_vec; // only used with vertex time uncertainties
106 #endif
107  // --- temporary numbers, used during update
108  std::vector<double> exp_arg_vec;
109  std::vector<double> exp_vec;
110  std::vector<double> sw_vec;
111  std::vector<double> swz_vec;
112  std::vector<double> swt_vec;
113  std::vector<double> se_vec;
114  std::vector<double> nuz_vec;
115  std::vector<double> nut_vec;
116  std::vector<double> szz_vec;
117  std::vector<double> stt_vec;
118  std::vector<double> szt_vec;
119 
120  unsigned int getSize() const { return zvtx_vec.size(); }
121 
122  void addItem(double new_zvtx, double new_tvtx, double new_rho) {
123  zvtx_vec.push_back(new_zvtx);
124  tvtx_vec.push_back(new_tvtx);
125  rho_vec.push_back(new_rho);
126  exp_arg_vec.push_back(0.0);
127  exp_vec.push_back(0.0);
128  swz_vec.push_back(0.0);
129  swt_vec.push_back(0.0);
130  se_vec.push_back(0.0);
131  nuz_vec.push_back(0.0);
132  nut_vec.push_back(0.0);
133  szz_vec.push_back(0.0);
134  stt_vec.push_back(0.0);
135  szt_vec.push_back(0.0);
136 #ifdef USEVTXDT2
137  dt2_vec.push_back(0.0);
138  sumw_vec.push_back(0.0);
139 #endif
140 
141  extractRaw();
142  }
143 
144  void insertItem(unsigned int k, double new_zvtx, double new_tvtx, double new_rho, track_t &tks) {
145  zvtx_vec.insert(zvtx_vec.begin() + k, new_zvtx);
146  tvtx_vec.insert(tvtx_vec.begin() + k, new_tvtx);
147  rho_vec.insert(rho_vec.begin() + k, new_rho);
148  exp_arg_vec.insert(exp_arg_vec.begin() + k, 0.0);
149  exp_vec.insert(exp_vec.begin() + k, 0.0);
150  swz_vec.insert(swz_vec.begin() + k, 0.0);
151  swt_vec.insert(swt_vec.begin() + k, 0.0);
152  se_vec.insert(se_vec.begin() + k, 0.0);
153  nuz_vec.insert(nuz_vec.begin() + k, 0.0);
154  nut_vec.insert(nut_vec.begin() + k, 0.0);
155  szz_vec.insert(szz_vec.begin() + k, 0.0);
156  stt_vec.insert(stt_vec.begin() + k, 0.0);
157  szt_vec.insert(szt_vec.begin() + k, 0.0);
158 #ifdef USEVTXDT2
159  dt2_vec.insert(dt2_vec.begin() + k, 0.0);
160  sumw_vec.insert(sumw_vec.begin() + k, 0.0);
161 #endif
162 
163  // adjust vertex lists of tracks
164  for (unsigned int i = 0; i < tks.getSize(); i++) {
165  if (tks.kmin[i] > k) {
166  tks.kmin[i]++;
167  }
168  if ((tks.kmax[i] >= k) || (tks.kmax[i] == tks.kmin[i])) {
169  tks.kmax[i]++;
170  }
171  }
172 
173  extractRaw();
174  }
175 
176  void removeItem(unsigned int k, track_t &tks) {
177  zvtx_vec.erase(zvtx_vec.begin() + k);
178  tvtx_vec.erase(tvtx_vec.begin() + k);
179  rho_vec.erase(rho_vec.begin() + k);
180  exp_arg_vec.erase(exp_arg_vec.begin() + k);
181  exp_vec.erase(exp_vec.begin() + k);
182  swz_vec.erase(swz_vec.begin() + k);
183  swt_vec.erase(swt_vec.begin() + k);
184  se_vec.erase(se_vec.begin() + k);
185  nuz_vec.erase(nuz_vec.begin() + k);
186  nut_vec.erase(nut_vec.begin() + k);
187  szz_vec.erase(szz_vec.begin() + k);
188  stt_vec.erase(stt_vec.begin() + k);
189  szt_vec.erase(szt_vec.begin() + k);
190 #ifdef USEVTXDT2
191  dt2_vec.erase(dt2_vec.begin() + k);
192  sumw_vec.erase(sumw_vec.begin() + k);
193 #endif
194 
195  // adjust vertex lists of tracks
196  for (unsigned int i = 0; i < tks.getSize(); i++) {
197  if (tks.kmax[i] > k) {
198  tks.kmax[i]--;
199  }
200  if ((tks.kmin[i] > k) || (((tks.kmax[i] < (tks.kmin[i] + 1)) && (tks.kmin[i] > 0)))) {
201  tks.kmin[i]--;
202  }
203  }
204 
205  extractRaw();
206  }
207 
208  unsigned int insertOrdered(double zvtx, double tvtx, double rho, track_t &tks) {
209  // insert a new cluster according to it's z-position, return the index at which it was inserted
210 
211  unsigned int k = 0;
212  for (; k < getSize(); k++) {
213  if (zvtx < zvtx_vec[k])
214  break;
215  }
216  insertItem(k, zvtx, tvtx, rho, tks);
217  return k;
218  }
219 
220  // pointers to the first element of vectors, needed for vectorized code
221  double *__restrict__ zvtx;
222  double *__restrict__ tvtx;
223  double *__restrict__ rho;
224  double *__restrict__ exp_arg;
225  double *__restrict__ exp;
226  double *__restrict__ swt;
227  double *__restrict__ swz;
228  double *__restrict__ se;
229  double *__restrict__ nuz;
230  double *__restrict__ nut;
231  double *__restrict__ szz;
232  double *__restrict__ stt;
233  double *__restrict__ szt;
234 #ifdef USEVTXDT2
235  double *__restrict__ dt2;
236  double *__restrict__ sumw;
237 #endif
238 
239  // has to be called everytime the items are modified
240  void extractRaw() {
241  zvtx = &zvtx_vec.front();
242  tvtx = &tvtx_vec.front();
243  rho = &rho_vec.front();
244  exp_arg = &exp_arg_vec.front();
245  exp = &exp_vec.front();
246  swz = &swz_vec.front();
247  swt = &swt_vec.front();
248  se = &se_vec.front();
249  nuz = &nuz_vec.front();
250  nut = &nut_vec.front();
251  szz = &szz_vec.front();
252  stt = &stt_vec.front();
253  szt = &szt_vec.front();
254 #ifdef USEVTXDT2
255  dt2 = &dt2_vec.front();
256  sumw = &sumw_vec.front();
257 #endif
258  }
259  };
260 
262 
264 
265  std::vector<std::vector<reco::TransientTrack> > clusterize(
266  const std::vector<reco::TransientTrack> &tracks) const override;
267 
268  std::vector<TransientVertex> vertices(const std::vector<reco::TransientTrack> &tracks) const override;
269 
270  track_t fill(const std::vector<reco::TransientTrack> &tracks) const;
271 
272  void set_vtx_range(double beta, track_t &gtracks, vertex_t &gvertices) const;
273 
274  void clear_vtx_range(track_t &gtracks, vertex_t &gvertices) const;
275 
276  unsigned int thermalize(
277  double beta, track_t &gtracks, vertex_t &gvertices, const double delta_max, const double rho0 = 0.) const;
278 
279  double update(
280  double beta, track_t &gtracks, vertex_t &gvertices, const double rho0 = 0, const bool updateTc = false) const;
281 
282  void dump(const double beta, const vertex_t &y, const track_t &tks, const int verbosity = 0) const;
283  bool zorder(vertex_t &y) const;
284  bool find_nearest(double z, double t, vertex_t &y, unsigned int &k_min, double dz, double dt) const;
285  bool merge(vertex_t &, track_t &, double &beta) const;
286  bool purge(vertex_t &, track_t &, double &, const double) const;
287  bool split(const double beta, track_t &t, vertex_t &y, double threshold = 1.) const;
288 
289  double beta0(const double betamax, track_t const &tks, vertex_t const &y) const;
290 
291  double get_Tc(const vertex_t &y, int k) const;
292  void verify(const vertex_t &v, const track_t &tks, unsigned int nv = 999999, unsigned int nt = 999999) const;
293 
294 private:
295  double zdumpcenter_;
296  double zdumpwidth_;
297 
298  double vertexSize_;
300  unsigned int maxIterations_;
302  double betamax_;
303  double betastop_;
304  double dzCutOff_;
305  double d0CutOff_;
306  double dtCutOff_;
307  double t0Max_;
308 
312  double zmerge_;
313  double tmerge_;
314  double betapurge_;
315 
316  unsigned int convergence_mode_;
317  double delta_highT_;
318  double delta_lowT_;
319 
320  double sel_zrange_;
321  const double zrange_min_ = 0.1; // smallest z-range to be included in a tracks cluster list
322 };
323 
324 //#ifndef DAClusterizerInZT_vect_h
325 #endif
float dt
Definition: AMPTWrapper.h:136
void removeItem(unsigned int k, track_t &tks)
std::vector< unsigned int > kmin
double get_Tc(const vertex_t &y, int k) const
DAClusterizerInZT_vect(const edm::ParameterSet &conf)
bool purge(vertex_t &, track_t &, double &, const double) const
void addItem(double new_zpca, double new_tpca, double new_dz2, double new_dt2, const reco::TransientTrack *new_tt, double new_tkwt)
std::vector< unsigned int > kmax
double update(double beta, track_t &gtracks, vertex_t &gvertices, const double rho0=0, const bool updateTc=false) const
void insertItem(unsigned int k, double new_zvtx, double new_tvtx, double new_rho, track_t &tks)
Definition: TTTypes.h:54
bool split(const double beta, track_t &t, vertex_t &y, double threshold=1.) const
double beta0(const double betamax, track_t const &tks, vertex_t const &y) const
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const override
void set_vtx_range(double beta, track_t &gtracks, vertex_t &gvertices) const
void addItem(double new_zvtx, double new_tvtx, double new_rho)
bool zorder(vertex_t &y) const
bool merge(vertex_t &, track_t &, double &beta) const
track_t fill(const std::vector< reco::TransientTrack > &tracks) const
void verify(const vertex_t &v, const track_t &tks, unsigned int nv=999999, unsigned int nt=999999) const
int nt
Definition: AMPTWrapper.h:42
const int verbosity
void insertItem(unsigned int i, double new_zpca, double new_tpca, double new_dz2, double new_dt2, const reco::TransientTrack *new_tt, double new_tkwt)
static void fillPSetDescription(edm::ParameterSetDescription &desc)
void dump(const double beta, const vertex_t &y, const track_t &tks, const int verbosity=0) const
std::vector< const reco::TransientTrack * > tt
std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const override
unsigned int thermalize(double beta, track_t &gtracks, vertex_t &gvertices, const double delta_max, const double rho0=0.) const
bool find_nearest(double z, double t, vertex_t &y, unsigned int &k_min, double dz, double dt) const
unsigned int insertOrdered(double zvtx, double tvtx, double rho, track_t &tks)
void clear_vtx_range(track_t &gtracks, vertex_t &gvertices) const