CMS 3D CMS Logo

HemisphereAlgo.cc
Go to the documentation of this file.
2 
4 
7 
8 using namespace std;
9 using std::vector;
10 
11 HemisphereAlgo::HemisphereAlgo(const std::vector<reco::CandidatePtr>& componentPtr_,
12  const int seed_method,
13  const int hemisphere_association_method)
14  : Object(componentPtr_),
15  Object_Group(),
16  Axis1(),
17  Axis2(),
18  seed_meth(seed_method),
19  hemi_meth(hemisphere_association_method),
20  status(0) {
21  for (int i = 0; i < (int)Object.size(); i++) {
22  Object_Noseed.push_back(0);
23  }
24 }
25 
26 vector<float> HemisphereAlgo::getAxis1() {
27  if (status != 1) {
28  this->reconstruct();
29  }
30  return Axis1;
31 }
32 vector<float> HemisphereAlgo::getAxis2() {
33  if (status != 1) {
34  this->reconstruct();
35  }
36  return Axis2;
37 }
38 
40  if (status != 1) {
41  this->reconstruct();
42  }
43  return Object_Group;
44 }
45 
47  int vsize = (int)Object.size();
48 
49  LogDebug("HemisphereAlgo") << " HemisphereAlgo method ";
50 
51  Object_Group.clear();
52  Axis1.clear();
53  Axis2.clear();
54 
55  for (int j = 0; j < vsize; j++) {
56  Object_Group.push_back(0);
57  }
58 
59  for (int j = 0; j < 5; j++) {
60  Axis1.push_back(0);
61  Axis2.push_back(0);
62  }
63 
64  for (int i = 0; i < vsize; i++) {
65  if ((*(Object)[i]).p() > (*Object[i]).energy() + 0.001) {
66  edm::LogWarning("HemisphereAlgo") << "Object " << i << " has E = " << (*Object[i]).energy()
67  << " less than P = " << (*Object[i]).p();
68  }
69  }
70 
71  LogDebug("HemisphereAlgo") << " Seeding method = " << seed_meth;
72  int I_Max = -1;
73  int J_Max = -1;
74 
75  if (seed_meth == 1) {
76  float P_Max = 0.;
77  float DeltaRP_Max = 0.;
78 
79  // take highest momentum object as first axis
80  for (int i = 0; i < vsize; i++) {
81  Object_Group[i] = 0;
82  if (Object_Noseed[i] == 0 && P_Max < (*Object[i]).p()) {
83  P_Max = (*Object[i]).p();
84  I_Max = i;
85  }
86  }
87 
88  Axis1[0] = (*Object[I_Max]).px() / (*Object[I_Max]).p();
89  Axis1[1] = (*Object[I_Max]).py() / (*Object[I_Max]).p();
90  Axis1[2] = (*Object[I_Max]).pz() / (*Object[I_Max]).p();
91  Axis1[3] = (*Object[I_Max]).p();
92  Axis1[4] = (*Object[I_Max]).energy();
93 
94  // take as second axis
95  for (int i = 0; i < vsize; i++) {
96  float DeltaR = deltaR((*Object[i]).eta(), (*Object[i]).phi(), (*Object[I_Max]).eta(), (*Object[I_Max]).phi());
97 
98  if (Object_Noseed[i] == 0 && DeltaR > 0.5) {
99  float DeltaRP = DeltaR * (*Object[i]).p();
100  if (DeltaRP > DeltaRP_Max) {
101  DeltaRP_Max = DeltaRP;
102  J_Max = i;
103  }
104  }
105  }
106 
107  if (J_Max >= 0) {
108  Axis2[0] = (*Object[J_Max]).px() / (*Object[J_Max]).p();
109  Axis2[1] = (*Object[J_Max]).py() / (*Object[J_Max]).p();
110  Axis2[2] = (*Object[J_Max]).pz() / (*Object[J_Max]).p();
111  Axis2[3] = (*Object[J_Max]).p();
112  Axis2[4] = (*Object[J_Max]).energy();
113 
114  } else {
115  return 0;
116  }
117  LogDebug("HemisphereAlgo") << " Axis 1 is Object = " << I_Max;
118  LogDebug("HemisphereAlgo") << " Axis 2 is Object = " << J_Max;
119 
120  } else if ((seed_meth == 2) | (seed_meth == 3)) {
121  float Mass_Max = 0.;
122  float InvariantMass = 0.;
123 
124  // maximize the invariant mass of two objects
125  for (int i = 0; i < vsize; i++) {
126  Object_Group[i] = 0;
127  if (Object_Noseed[i] == 0) {
128  for (int j = i + 1; j < vsize; j++) {
129  if (Object_Noseed[j] == 0) {
130  // either the invariant mass
131  if (seed_meth == 2) {
132  InvariantMass =
133  ((*Object[i]).energy() + (*Object[j]).energy()) * ((*Object[i]).energy() + (*Object[j]).energy()) -
134  ((*Object[i]).px() + (*Object[j]).px()) * ((*Object[i]).px() + (*Object[j]).px()) -
135  ((*Object[i]).py() + (*Object[j]).py()) * ((*Object[i]).py() + (*Object[j]).py()) -
136  ((*Object[i]).pz() + (*Object[j]).pz()) * ((*Object[i]).pz() + (*Object[j]).pz());
137  }
138  // or the transverse mass
139  else if (seed_meth == 3) {
140  float pti = sqrt((*Object[i]).px() * (*Object[i]).px() + (*Object[i]).py() * (*Object[i]).py());
141  float ptj = sqrt((*Object[j]).px() * (*Object[j]).px() + (*Object[j]).py() * (*Object[j]).py());
142  InvariantMass =
143  2. * (pti * ptj - (*Object[i]).px() * (*Object[j]).px() - (*Object[i]).py() * (*Object[j]).py());
144  }
145  if (Mass_Max < InvariantMass) {
146  Mass_Max = InvariantMass;
147  I_Max = i;
148  J_Max = j;
149  }
150  }
151  }
152  }
153  }
154 
155  if (J_Max > 0) {
156  Axis1[0] = (*Object[I_Max]).px() / (*Object[I_Max]).p();
157  Axis1[1] = (*Object[I_Max]).py() / (*Object[I_Max]).p();
158  Axis1[2] = (*Object[I_Max]).pz() / (*Object[I_Max]).p();
159 
160  Axis1[3] = (*Object[I_Max]).p();
161  Axis1[4] = (*Object[I_Max]).energy();
162 
163  Axis2[0] = (*Object[J_Max]).px() / (*Object[J_Max]).p();
164  Axis2[1] = (*Object[J_Max]).py() / (*Object[J_Max]).p();
165  Axis2[2] = (*Object[J_Max]).pz() / (*Object[J_Max]).p();
166 
167  Axis2[3] = (*Object[J_Max]).p();
168  Axis2[4] = (*Object[J_Max]).energy();
169 
170  } else {
171  return 0;
172  }
173  LogDebug("HemisphereAlgo") << " Axis 1 is Object = " << I_Max;
174  LogDebug("HemisphereAlgo") << " Axis 2 is Object = " << J_Max;
175 
176  } else {
177  throw cms::Exception("Configuration") << "Please give a valid seeding method!";
178  }
179 
180  // seeding done
181  // now do the hemisphere association
182 
183  LogDebug("HemisphereAlgo") << " Association method = " << hemi_meth;
184 
185  int numLoop = 0;
186  bool I_Move = true;
187 
188  while (I_Move && (numLoop < 100)) {
189  I_Move = false;
190  numLoop++;
191  LogDebug("HemisphereAlgo") << " Iteration = " << numLoop;
192 
193  float Sum1_Px = 0.;
194  float Sum1_Py = 0.;
195  float Sum1_Pz = 0.;
196  float Sum1_E = 0.;
197  float Sum2_Px = 0.;
198  float Sum2_Py = 0.;
199  float Sum2_Pz = 0.;
200  float Sum2_E = 0.;
201 
202  if (hemi_meth == 1) {
203  for (int i = 0; i < vsize; i++) {
204  float P_Long1 = (*Object[i]).px() * Axis1[0] + (*Object[i]).py() * Axis1[1] + (*Object[i]).pz() * Axis1[2];
205  float P_Long2 = (*Object[i]).px() * Axis2[0] + (*Object[i]).py() * Axis2[1] + (*Object[i]).pz() * Axis2[2];
206  if (P_Long1 >= P_Long2) {
207  if (Object_Group[i] != 1) {
208  I_Move = true;
209  }
210  Object_Group[i] = 1;
211  Sum1_Px += (*Object[i]).px();
212  Sum1_Py += (*Object[i]).py();
213  Sum1_Pz += (*Object[i]).pz();
214  Sum1_E += (*Object[i]).energy();
215  } else {
216  if (Object_Group[i] != 2) {
217  I_Move = true;
218  }
219  Object_Group[i] = 2;
220  Sum2_Px += (*Object[i]).px();
221  Sum2_Py += (*Object[i]).py();
222  Sum2_Pz += (*Object[i]).pz();
223  Sum2_E += (*Object[i]).energy();
224  }
225  }
226 
227  } else if (hemi_meth == 2 || hemi_meth == 3) {
228  for (int i = 0; i < vsize; i++) {
229  if (i == I_Max) {
230  Object_Group[i] = 1;
231  Sum1_Px += (*Object[i]).px();
232  Sum1_Py += (*Object[i]).py();
233  Sum1_Pz += (*Object[i]).pz();
234  Sum1_E += (*Object[i]).energy();
235  } else if (i == J_Max) {
236  Object_Group[i] = 2;
237  Sum2_Px += (*Object[i]).px();
238  Sum2_Py += (*Object[i]).py();
239  Sum2_Pz += (*Object[i]).pz();
240  Sum2_E += (*Object[i]).energy();
241  } else {
242  if (!I_Move) {
243  float NewAxis1_Px = Axis1[0] * Axis1[3];
244  float NewAxis1_Py = Axis1[1] * Axis1[3];
245  float NewAxis1_Pz = Axis1[2] * Axis1[3];
246  float NewAxis1_E = Axis1[4];
247 
248  float NewAxis2_Px = Axis2[0] * Axis2[3];
249  float NewAxis2_Py = Axis2[1] * Axis2[3];
250  float NewAxis2_Pz = Axis2[2] * Axis2[3];
251  float NewAxis2_E = Axis2[4];
252 
253  if (Object_Group[i] == 1) {
254  NewAxis1_Px = NewAxis1_Px - (*Object[i]).px();
255  NewAxis1_Py = NewAxis1_Py - (*Object[i]).py();
256  NewAxis1_Pz = NewAxis1_Pz - (*Object[i]).pz();
257  NewAxis1_E = NewAxis1_E - (*Object[i]).energy();
258 
259  } else if (Object_Group[i] == 2) {
260  NewAxis2_Px = NewAxis2_Px - (*Object[i]).px();
261  NewAxis2_Py = NewAxis2_Py - (*Object[i]).py();
262  NewAxis2_Pz = NewAxis2_Pz - (*Object[i]).pz();
263  NewAxis2_E = NewAxis2_E - (*Object[i]).energy();
264  }
265 
266  float mass1 =
267  NewAxis1_E -
268  (((*Object[i]).px() * NewAxis1_Px + (*Object[i]).py() * NewAxis1_Py + (*Object[i]).pz() * NewAxis1_Pz) /
269  (*Object[i]).p());
270 
271  float mass2 =
272  NewAxis2_E -
273  (((*Object[i]).px() * NewAxis2_Px + (*Object[i]).py() * NewAxis2_Py + (*Object[i]).pz() * NewAxis2_Pz) /
274  (*Object[i]).p());
275 
276  if (hemi_meth == 3) {
277  mass1 *= NewAxis1_E / ((NewAxis1_E + (*Object[i]).energy()) * (NewAxis1_E + (*Object[i]).energy()));
278 
279  mass2 *= NewAxis2_E / ((NewAxis2_E + (*Object[i]).energy()) * (NewAxis2_E + (*Object[i]).energy()));
280  }
281 
282  if (mass1 < mass2) {
283  if (Object_Group[i] != 1) {
284  I_Move = true;
285  }
286  Object_Group[i] = 1;
287 
288  Sum1_Px += (*Object[i]).px();
289  Sum1_Py += (*Object[i]).py();
290  Sum1_Pz += (*Object[i]).pz();
291  Sum1_E += (*Object[i]).energy();
292  } else {
293  if (Object_Group[i] != 2) {
294  I_Move = true;
295  }
296  Object_Group[i] = 2;
297  Sum2_Px += (*Object[i]).px();
298  Sum2_Py += (*Object[i]).py();
299  Sum2_Pz += (*Object[i]).pz();
300  Sum2_E += (*Object[i]).energy();
301  }
302 
303  } else {
304  if (Object_Group[i] == 1) {
305  Sum1_Px += (*Object[i]).px();
306  Sum1_Py += (*Object[i]).py();
307  Sum1_Pz += (*Object[i]).pz();
308  Sum1_E += (*Object[i]).energy();
309  }
310  if (Object_Group[i] == 2) {
311  Sum2_Px += (*Object[i]).px();
312  Sum2_Py += (*Object[i]).py();
313  Sum2_Pz += (*Object[i]).pz();
314  Sum2_E += (*Object[i]).energy();
315  }
316  }
317  }
318  }
319  } else {
320  throw cms::Exception("Configuration") << "Please give a valid hemisphere association method!";
321  }
322 
323  // recomputing the axes
324 
325  Axis1[3] = sqrt(Sum1_Px * Sum1_Px + Sum1_Py * Sum1_Py + Sum1_Pz * Sum1_Pz);
326  if (Axis1[3] < 0.0001) {
327  edm::LogWarning("HemisphereAlgo") << "ZERO objects in group 1! ";
328  } else {
329  Axis1[0] = Sum1_Px / Axis1[3];
330  Axis1[1] = Sum1_Py / Axis1[3];
331  Axis1[2] = Sum1_Pz / Axis1[3];
332  Axis1[4] = Sum1_E;
333  }
334 
335  Axis2[3] = sqrt(Sum2_Px * Sum2_Px + Sum2_Py * Sum2_Py + Sum2_Pz * Sum2_Pz);
336  if (Axis2[3] < 0.0001) {
337  edm::LogWarning("HemisphereAlgo") << "ZERO objects in group 2!";
338  } else {
339  Axis2[0] = Sum2_Px / Axis2[3];
340  Axis2[1] = Sum2_Py / Axis2[3];
341  Axis2[2] = Sum2_Pz / Axis2[3];
342  Axis2[4] = Sum2_E;
343  }
344 
345  LogDebug("HemisphereAlgo") << " Grouping = ";
346  for (int i = 0; i < vsize; i++) {
347  LogTrace("HemisphereAlgo") << " " << Object_Group[i];
348  }
349  LogTrace("HemisphereAlgo") << std::endl;
350  }
351 
352  status = 1;
353  return 1;
354 }
Definition: DeltaR.py:1
std::vector< int > Object_Noseed
std::vector< float > Axis2
std::vector< float > Axis1
std::vector< int > Object_Group
std::vector< float > getAxis1()
#define LogTrace(id)
HemisphereAlgo(const std::vector< reco::CandidatePtr > &componentRefs_, const int seed_method=0, const int hemisphere_association_method=0)
T sqrt(T t)
Definition: SSEVec.h:23
std::vector< int > getGrouping()
std::vector< reco::CandidatePtr > Object
Log< level::Warning, false > LogWarning
std::vector< float > getAxis2()
#define LogDebug(id)