SFS: atlas_susy_2018_31.cpp

File atlas_susy_2018_31.cpp, 35.4 KB (added by Benjamin Fuks, 4 years ago)
Line 
1#include "SampleAnalyzer/User/Analyzer/atlas_susy_2018_31.h"
2// #include "SampleAnalyzer/Commons/Vector/MAMatrix.h"
3using namespace MA5;
4using namespace std;
5
6template<typename T1, typename T2> std::vector<const T1*>
7 Removal(std::vector<const T1*>&, std::vector<const T2*>&, const double &);
8
9// static void rotateXY(MA5::MAMatrix &mat, MA5::MAMatrix &mat_new,MAdouble64 phi);
10// static void add2x2Matrices(MAMatrix &mat1, MAMatrix &mat2);
11// -----------------------------------------------------------------------------
12// Initialize
13// function called one time at the beginning of the analysis
14// -----------------------------------------------------------------------------
15bool atlas_susy_2018_31::Initialize(const MA5::Configuration& cfg,
16 const std::map<std::string,std::string>& parameters)
17{
18 INFO << " <><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><>" << endmsg;
19 INFO << " <> ATLAS SUSY 2018 031 <>" << endmsg;
20 INFO << " <> sbottom, multi-Bjets + MET @ 13 TeV, 139 fb^-1 <>" << endmsg;
21 INFO << " <> arXiv:1908.03122 <>" << endmsg;
22 INFO << " <> Recasted by : J.Y. Araz, B. Fuks <>" << endmsg;
23 INFO << " <> Contact : jack.araz@durham.ac.uk,fuks@lpthe.jussieu.fr <>" << endmsg;
24 INFO << " <> Based on MadAnalysis 5 v1.9.4 and above <>" << endmsg;
25 INFO << " <> DOI: XX.YYYY/ZZZ <>" << endmsg;
26 INFO << " <> For more information, see <>" << endmsg;
27 INFO << " <> http://madanalysis.irmp.ucl.ac.be/wiki/SFS <>" << endmsg;
28 INFO << " <><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><><>" << endmsg;
29
30 // Initializing PhysicsService for MC
31 PHYSICS->mcConfig().Reset();
32
33 // definition of the multiparticle "hadronic"
34 PHYSICS->mcConfig().AddHadronicId(-20543);
35 PHYSICS->mcConfig().AddHadronicId(-20533);
36 PHYSICS->mcConfig().AddHadronicId(-20523);
37 PHYSICS->mcConfig().AddHadronicId(-20513);
38 PHYSICS->mcConfig().AddHadronicId(-20433);
39 PHYSICS->mcConfig().AddHadronicId(-20423);
40 PHYSICS->mcConfig().AddHadronicId(-20413);
41 PHYSICS->mcConfig().AddHadronicId(-20323);
42 PHYSICS->mcConfig().AddHadronicId(-20313);
43 PHYSICS->mcConfig().AddHadronicId(-20213);
44 PHYSICS->mcConfig().AddHadronicId(-10543);
45 PHYSICS->mcConfig().AddHadronicId(-10541);
46 PHYSICS->mcConfig().AddHadronicId(-10533);
47 PHYSICS->mcConfig().AddHadronicId(-10531);
48 PHYSICS->mcConfig().AddHadronicId(-10523);
49 PHYSICS->mcConfig().AddHadronicId(-10521);
50 PHYSICS->mcConfig().AddHadronicId(-10513);
51 PHYSICS->mcConfig().AddHadronicId(-10511);
52 PHYSICS->mcConfig().AddHadronicId(-10433);
53 PHYSICS->mcConfig().AddHadronicId(-10431);
54 PHYSICS->mcConfig().AddHadronicId(-10423);
55 PHYSICS->mcConfig().AddHadronicId(-10421);
56 PHYSICS->mcConfig().AddHadronicId(-10413);
57 PHYSICS->mcConfig().AddHadronicId(-10411);
58 PHYSICS->mcConfig().AddHadronicId(-10323);
59 PHYSICS->mcConfig().AddHadronicId(-10321);
60 PHYSICS->mcConfig().AddHadronicId(-10313);
61 PHYSICS->mcConfig().AddHadronicId(-10311);
62 PHYSICS->mcConfig().AddHadronicId(-10213);
63 PHYSICS->mcConfig().AddHadronicId(-10211);
64 PHYSICS->mcConfig().AddHadronicId(-5554);
65 PHYSICS->mcConfig().AddHadronicId(-5544);
66 PHYSICS->mcConfig().AddHadronicId(-5542);
67 PHYSICS->mcConfig().AddHadronicId(-5534);
68 PHYSICS->mcConfig().AddHadronicId(-5532);
69 PHYSICS->mcConfig().AddHadronicId(-5524);
70 PHYSICS->mcConfig().AddHadronicId(-5522);
71 PHYSICS->mcConfig().AddHadronicId(-5514);
72 PHYSICS->mcConfig().AddHadronicId(-5512);
73 PHYSICS->mcConfig().AddHadronicId(-5503);
74 PHYSICS->mcConfig().AddHadronicId(-5444);
75 PHYSICS->mcConfig().AddHadronicId(-5442);
76 PHYSICS->mcConfig().AddHadronicId(-5434);
77 PHYSICS->mcConfig().AddHadronicId(-5432);
78 PHYSICS->mcConfig().AddHadronicId(-5424);
79 PHYSICS->mcConfig().AddHadronicId(-5422);
80 PHYSICS->mcConfig().AddHadronicId(-5414);
81 PHYSICS->mcConfig().AddHadronicId(-5412);
82 PHYSICS->mcConfig().AddHadronicId(-5403);
83 PHYSICS->mcConfig().AddHadronicId(-5401);
84 PHYSICS->mcConfig().AddHadronicId(-5342);
85 PHYSICS->mcConfig().AddHadronicId(-5334);
86 PHYSICS->mcConfig().AddHadronicId(-5332);
87 PHYSICS->mcConfig().AddHadronicId(-5324);
88 PHYSICS->mcConfig().AddHadronicId(-5322);
89 PHYSICS->mcConfig().AddHadronicId(-5314);
90 PHYSICS->mcConfig().AddHadronicId(-5312);
91 PHYSICS->mcConfig().AddHadronicId(-5303);
92 PHYSICS->mcConfig().AddHadronicId(-5301);
93 PHYSICS->mcConfig().AddHadronicId(-5242);
94 PHYSICS->mcConfig().AddHadronicId(-5232);
95 PHYSICS->mcConfig().AddHadronicId(-5224);
96 PHYSICS->mcConfig().AddHadronicId(-5222);
97 PHYSICS->mcConfig().AddHadronicId(-5214);
98 PHYSICS->mcConfig().AddHadronicId(-5212);
99 PHYSICS->mcConfig().AddHadronicId(-5203);
100 PHYSICS->mcConfig().AddHadronicId(-5201);
101 PHYSICS->mcConfig().AddHadronicId(-5142);
102 PHYSICS->mcConfig().AddHadronicId(-5132);
103 PHYSICS->mcConfig().AddHadronicId(-5122);
104 PHYSICS->mcConfig().AddHadronicId(-5114);
105 PHYSICS->mcConfig().AddHadronicId(-5112);
106 PHYSICS->mcConfig().AddHadronicId(-5103);
107 PHYSICS->mcConfig().AddHadronicId(-5101);
108 PHYSICS->mcConfig().AddHadronicId(-4444);
109 PHYSICS->mcConfig().AddHadronicId(-4434);
110 PHYSICS->mcConfig().AddHadronicId(-4432);
111 PHYSICS->mcConfig().AddHadronicId(-4424);
112 PHYSICS->mcConfig().AddHadronicId(-4422);
113 PHYSICS->mcConfig().AddHadronicId(-4414);
114 PHYSICS->mcConfig().AddHadronicId(-4412);
115 PHYSICS->mcConfig().AddHadronicId(-4403);
116 PHYSICS->mcConfig().AddHadronicId(-4334);
117 PHYSICS->mcConfig().AddHadronicId(-4332);
118 PHYSICS->mcConfig().AddHadronicId(-4324);
119 PHYSICS->mcConfig().AddHadronicId(-4322);
120 PHYSICS->mcConfig().AddHadronicId(-4314);
121 PHYSICS->mcConfig().AddHadronicId(-4312);
122 PHYSICS->mcConfig().AddHadronicId(-4303);
123 PHYSICS->mcConfig().AddHadronicId(-4301);
124 PHYSICS->mcConfig().AddHadronicId(-4232);
125 PHYSICS->mcConfig().AddHadronicId(-4224);
126 PHYSICS->mcConfig().AddHadronicId(-4222);
127 PHYSICS->mcConfig().AddHadronicId(-4214);
128 PHYSICS->mcConfig().AddHadronicId(-4212);
129 PHYSICS->mcConfig().AddHadronicId(-4203);
130 PHYSICS->mcConfig().AddHadronicId(-4201);
131 PHYSICS->mcConfig().AddHadronicId(-4132);
132 PHYSICS->mcConfig().AddHadronicId(-4122);
133 PHYSICS->mcConfig().AddHadronicId(-4114);
134 PHYSICS->mcConfig().AddHadronicId(-4112);
135 PHYSICS->mcConfig().AddHadronicId(-4103);
136 PHYSICS->mcConfig().AddHadronicId(-4101);
137 PHYSICS->mcConfig().AddHadronicId(-3334);
138 PHYSICS->mcConfig().AddHadronicId(-3324);
139 PHYSICS->mcConfig().AddHadronicId(-3322);
140 PHYSICS->mcConfig().AddHadronicId(-3314);
141 PHYSICS->mcConfig().AddHadronicId(-3312);
142 PHYSICS->mcConfig().AddHadronicId(-3303);
143 PHYSICS->mcConfig().AddHadronicId(-3224);
144 PHYSICS->mcConfig().AddHadronicId(-3222);
145 PHYSICS->mcConfig().AddHadronicId(-3214);
146 PHYSICS->mcConfig().AddHadronicId(-3212);
147 PHYSICS->mcConfig().AddHadronicId(-3203);
148 PHYSICS->mcConfig().AddHadronicId(-3201);
149 PHYSICS->mcConfig().AddHadronicId(-3122);
150 PHYSICS->mcConfig().AddHadronicId(-3114);
151 PHYSICS->mcConfig().AddHadronicId(-3112);
152 PHYSICS->mcConfig().AddHadronicId(-3103);
153 PHYSICS->mcConfig().AddHadronicId(-3101);
154 PHYSICS->mcConfig().AddHadronicId(-2224);
155 PHYSICS->mcConfig().AddHadronicId(-2214);
156 PHYSICS->mcConfig().AddHadronicId(-2212);
157 PHYSICS->mcConfig().AddHadronicId(-2203);
158 PHYSICS->mcConfig().AddHadronicId(-2114);
159 PHYSICS->mcConfig().AddHadronicId(-2112);
160 PHYSICS->mcConfig().AddHadronicId(-2103);
161 PHYSICS->mcConfig().AddHadronicId(-2101);
162 PHYSICS->mcConfig().AddHadronicId(-1114);
163 PHYSICS->mcConfig().AddHadronicId(-1103);
164 PHYSICS->mcConfig().AddHadronicId(-545);
165 PHYSICS->mcConfig().AddHadronicId(-543);
166 PHYSICS->mcConfig().AddHadronicId(-541);
167 PHYSICS->mcConfig().AddHadronicId(-535);
168 PHYSICS->mcConfig().AddHadronicId(-533);
169 PHYSICS->mcConfig().AddHadronicId(-531);
170 PHYSICS->mcConfig().AddHadronicId(-525);
171 PHYSICS->mcConfig().AddHadronicId(-523);
172 PHYSICS->mcConfig().AddHadronicId(-521);
173 PHYSICS->mcConfig().AddHadronicId(-515);
174 PHYSICS->mcConfig().AddHadronicId(-513);
175 PHYSICS->mcConfig().AddHadronicId(-511);
176 PHYSICS->mcConfig().AddHadronicId(-435);
177 PHYSICS->mcConfig().AddHadronicId(-433);
178 PHYSICS->mcConfig().AddHadronicId(-431);
179 PHYSICS->mcConfig().AddHadronicId(-425);
180 PHYSICS->mcConfig().AddHadronicId(-423);
181 PHYSICS->mcConfig().AddHadronicId(-421);
182 PHYSICS->mcConfig().AddHadronicId(-415);
183 PHYSICS->mcConfig().AddHadronicId(-413);
184 PHYSICS->mcConfig().AddHadronicId(-411);
185 PHYSICS->mcConfig().AddHadronicId(-325);
186 PHYSICS->mcConfig().AddHadronicId(-323);
187 PHYSICS->mcConfig().AddHadronicId(-321);
188 PHYSICS->mcConfig().AddHadronicId(-315);
189 PHYSICS->mcConfig().AddHadronicId(-313);
190 PHYSICS->mcConfig().AddHadronicId(-311);
191 PHYSICS->mcConfig().AddHadronicId(-215);
192 PHYSICS->mcConfig().AddHadronicId(-213);
193 PHYSICS->mcConfig().AddHadronicId(-211);
194 PHYSICS->mcConfig().AddHadronicId(111);
195 PHYSICS->mcConfig().AddHadronicId(113);
196 PHYSICS->mcConfig().AddHadronicId(115);
197 PHYSICS->mcConfig().AddHadronicId(130);
198 PHYSICS->mcConfig().AddHadronicId(211);
199 PHYSICS->mcConfig().AddHadronicId(213);
200 PHYSICS->mcConfig().AddHadronicId(215);
201 PHYSICS->mcConfig().AddHadronicId(221);
202 PHYSICS->mcConfig().AddHadronicId(223);
203 PHYSICS->mcConfig().AddHadronicId(225);
204 PHYSICS->mcConfig().AddHadronicId(310);
205 PHYSICS->mcConfig().AddHadronicId(311);
206 PHYSICS->mcConfig().AddHadronicId(313);
207 PHYSICS->mcConfig().AddHadronicId(315);
208 PHYSICS->mcConfig().AddHadronicId(321);
209 PHYSICS->mcConfig().AddHadronicId(323);
210 PHYSICS->mcConfig().AddHadronicId(325);
211 PHYSICS->mcConfig().AddHadronicId(331);
212 PHYSICS->mcConfig().AddHadronicId(333);
213 PHYSICS->mcConfig().AddHadronicId(335);
214 PHYSICS->mcConfig().AddHadronicId(411);
215 PHYSICS->mcConfig().AddHadronicId(413);
216 PHYSICS->mcConfig().AddHadronicId(415);
217 PHYSICS->mcConfig().AddHadronicId(421);
218 PHYSICS->mcConfig().AddHadronicId(423);
219 PHYSICS->mcConfig().AddHadronicId(425);
220 PHYSICS->mcConfig().AddHadronicId(431);
221 PHYSICS->mcConfig().AddHadronicId(433);
222 PHYSICS->mcConfig().AddHadronicId(435);
223 PHYSICS->mcConfig().AddHadronicId(441);
224 PHYSICS->mcConfig().AddHadronicId(443);
225 PHYSICS->mcConfig().AddHadronicId(445);
226 PHYSICS->mcConfig().AddHadronicId(511);
227 PHYSICS->mcConfig().AddHadronicId(513);
228 PHYSICS->mcConfig().AddHadronicId(515);
229 PHYSICS->mcConfig().AddHadronicId(521);
230 PHYSICS->mcConfig().AddHadronicId(523);
231 PHYSICS->mcConfig().AddHadronicId(525);
232 PHYSICS->mcConfig().AddHadronicId(531);
233 PHYSICS->mcConfig().AddHadronicId(533);
234 PHYSICS->mcConfig().AddHadronicId(535);
235 PHYSICS->mcConfig().AddHadronicId(541);
236 PHYSICS->mcConfig().AddHadronicId(543);
237 PHYSICS->mcConfig().AddHadronicId(545);
238 PHYSICS->mcConfig().AddHadronicId(551);
239 PHYSICS->mcConfig().AddHadronicId(553);
240 PHYSICS->mcConfig().AddHadronicId(555);
241 PHYSICS->mcConfig().AddHadronicId(1103);
242 PHYSICS->mcConfig().AddHadronicId(1114);
243 PHYSICS->mcConfig().AddHadronicId(2101);
244 PHYSICS->mcConfig().AddHadronicId(2103);
245 PHYSICS->mcConfig().AddHadronicId(2112);
246 PHYSICS->mcConfig().AddHadronicId(2114);
247 PHYSICS->mcConfig().AddHadronicId(2203);
248 PHYSICS->mcConfig().AddHadronicId(2212);
249 PHYSICS->mcConfig().AddHadronicId(2214);
250 PHYSICS->mcConfig().AddHadronicId(2224);
251 PHYSICS->mcConfig().AddHadronicId(3101);
252 PHYSICS->mcConfig().AddHadronicId(3103);
253 PHYSICS->mcConfig().AddHadronicId(3112);
254 PHYSICS->mcConfig().AddHadronicId(3114);
255 PHYSICS->mcConfig().AddHadronicId(3122);
256 PHYSICS->mcConfig().AddHadronicId(3201);
257 PHYSICS->mcConfig().AddHadronicId(3203);
258 PHYSICS->mcConfig().AddHadronicId(3212);
259 PHYSICS->mcConfig().AddHadronicId(3214);
260 PHYSICS->mcConfig().AddHadronicId(3222);
261 PHYSICS->mcConfig().AddHadronicId(3224);
262 PHYSICS->mcConfig().AddHadronicId(3303);
263 PHYSICS->mcConfig().AddHadronicId(3312);
264 PHYSICS->mcConfig().AddHadronicId(3314);
265 PHYSICS->mcConfig().AddHadronicId(3322);
266 PHYSICS->mcConfig().AddHadronicId(3324);
267 PHYSICS->mcConfig().AddHadronicId(3334);
268 PHYSICS->mcConfig().AddHadronicId(4101);
269 PHYSICS->mcConfig().AddHadronicId(4103);
270 PHYSICS->mcConfig().AddHadronicId(4112);
271 PHYSICS->mcConfig().AddHadronicId(4114);
272 PHYSICS->mcConfig().AddHadronicId(4122);
273 PHYSICS->mcConfig().AddHadronicId(4132);
274 PHYSICS->mcConfig().AddHadronicId(4201);
275 PHYSICS->mcConfig().AddHadronicId(4203);
276 PHYSICS->mcConfig().AddHadronicId(4212);
277 PHYSICS->mcConfig().AddHadronicId(4214);
278 PHYSICS->mcConfig().AddHadronicId(4222);
279 PHYSICS->mcConfig().AddHadronicId(4224);
280 PHYSICS->mcConfig().AddHadronicId(4232);
281 PHYSICS->mcConfig().AddHadronicId(4301);
282 PHYSICS->mcConfig().AddHadronicId(4303);
283 PHYSICS->mcConfig().AddHadronicId(4312);
284 PHYSICS->mcConfig().AddHadronicId(4314);
285 PHYSICS->mcConfig().AddHadronicId(4322);
286 PHYSICS->mcConfig().AddHadronicId(4324);
287 PHYSICS->mcConfig().AddHadronicId(4332);
288 PHYSICS->mcConfig().AddHadronicId(4334);
289 PHYSICS->mcConfig().AddHadronicId(4403);
290 PHYSICS->mcConfig().AddHadronicId(4412);
291 PHYSICS->mcConfig().AddHadronicId(4414);
292 PHYSICS->mcConfig().AddHadronicId(4422);
293 PHYSICS->mcConfig().AddHadronicId(4424);
294 PHYSICS->mcConfig().AddHadronicId(4432);
295 PHYSICS->mcConfig().AddHadronicId(4434);
296 PHYSICS->mcConfig().AddHadronicId(4444);
297 PHYSICS->mcConfig().AddHadronicId(5101);
298 PHYSICS->mcConfig().AddHadronicId(5103);
299 PHYSICS->mcConfig().AddHadronicId(5112);
300 PHYSICS->mcConfig().AddHadronicId(5114);
301 PHYSICS->mcConfig().AddHadronicId(5122);
302 PHYSICS->mcConfig().AddHadronicId(5132);
303 PHYSICS->mcConfig().AddHadronicId(5142);
304 PHYSICS->mcConfig().AddHadronicId(5201);
305 PHYSICS->mcConfig().AddHadronicId(5203);
306 PHYSICS->mcConfig().AddHadronicId(5212);
307 PHYSICS->mcConfig().AddHadronicId(5214);
308 PHYSICS->mcConfig().AddHadronicId(5222);
309 PHYSICS->mcConfig().AddHadronicId(5224);
310 PHYSICS->mcConfig().AddHadronicId(5232);
311 PHYSICS->mcConfig().AddHadronicId(5242);
312 PHYSICS->mcConfig().AddHadronicId(5301);
313 PHYSICS->mcConfig().AddHadronicId(5303);
314 PHYSICS->mcConfig().AddHadronicId(5312);
315 PHYSICS->mcConfig().AddHadronicId(5314);
316 PHYSICS->mcConfig().AddHadronicId(5322);
317 PHYSICS->mcConfig().AddHadronicId(5324);
318 PHYSICS->mcConfig().AddHadronicId(5332);
319 PHYSICS->mcConfig().AddHadronicId(5334);
320 PHYSICS->mcConfig().AddHadronicId(5342);
321 PHYSICS->mcConfig().AddHadronicId(5401);
322 PHYSICS->mcConfig().AddHadronicId(5403);
323 PHYSICS->mcConfig().AddHadronicId(5412);
324 PHYSICS->mcConfig().AddHadronicId(5414);
325 PHYSICS->mcConfig().AddHadronicId(5422);
326 PHYSICS->mcConfig().AddHadronicId(5424);
327 PHYSICS->mcConfig().AddHadronicId(5432);
328 PHYSICS->mcConfig().AddHadronicId(5434);
329 PHYSICS->mcConfig().AddHadronicId(5442);
330 PHYSICS->mcConfig().AddHadronicId(5444);
331 PHYSICS->mcConfig().AddHadronicId(5503);
332 PHYSICS->mcConfig().AddHadronicId(5512);
333 PHYSICS->mcConfig().AddHadronicId(5514);
334 PHYSICS->mcConfig().AddHadronicId(5522);
335 PHYSICS->mcConfig().AddHadronicId(5524);
336 PHYSICS->mcConfig().AddHadronicId(5532);
337 PHYSICS->mcConfig().AddHadronicId(5534);
338 PHYSICS->mcConfig().AddHadronicId(5542);
339 PHYSICS->mcConfig().AddHadronicId(5544);
340 PHYSICS->mcConfig().AddHadronicId(5554);
341 PHYSICS->mcConfig().AddHadronicId(10111);
342 PHYSICS->mcConfig().AddHadronicId(10113);
343 PHYSICS->mcConfig().AddHadronicId(10211);
344 PHYSICS->mcConfig().AddHadronicId(10213);
345 PHYSICS->mcConfig().AddHadronicId(10221);
346 PHYSICS->mcConfig().AddHadronicId(10223);
347 PHYSICS->mcConfig().AddHadronicId(10311);
348 PHYSICS->mcConfig().AddHadronicId(10313);
349 PHYSICS->mcConfig().AddHadronicId(10321);
350 PHYSICS->mcConfig().AddHadronicId(10323);
351 PHYSICS->mcConfig().AddHadronicId(10331);
352 PHYSICS->mcConfig().AddHadronicId(10333);
353 PHYSICS->mcConfig().AddHadronicId(10411);
354 PHYSICS->mcConfig().AddHadronicId(10413);
355 PHYSICS->mcConfig().AddHadronicId(10421);
356 PHYSICS->mcConfig().AddHadronicId(10423);
357 PHYSICS->mcConfig().AddHadronicId(10431);
358 PHYSICS->mcConfig().AddHadronicId(10433);
359 PHYSICS->mcConfig().AddHadronicId(10441);
360 PHYSICS->mcConfig().AddHadronicId(10443);
361 PHYSICS->mcConfig().AddHadronicId(10511);
362 PHYSICS->mcConfig().AddHadronicId(10513);
363 PHYSICS->mcConfig().AddHadronicId(10521);
364 PHYSICS->mcConfig().AddHadronicId(10523);
365 PHYSICS->mcConfig().AddHadronicId(10531);
366 PHYSICS->mcConfig().AddHadronicId(10533);
367 PHYSICS->mcConfig().AddHadronicId(10541);
368 PHYSICS->mcConfig().AddHadronicId(10543);
369 PHYSICS->mcConfig().AddHadronicId(10551);
370 PHYSICS->mcConfig().AddHadronicId(10553);
371 PHYSICS->mcConfig().AddHadronicId(20113);
372 PHYSICS->mcConfig().AddHadronicId(20213);
373 PHYSICS->mcConfig().AddHadronicId(20223);
374 PHYSICS->mcConfig().AddHadronicId(20313);
375 PHYSICS->mcConfig().AddHadronicId(20323);
376 PHYSICS->mcConfig().AddHadronicId(20333);
377 PHYSICS->mcConfig().AddHadronicId(20413);
378 PHYSICS->mcConfig().AddHadronicId(20423);
379 PHYSICS->mcConfig().AddHadronicId(20433);
380 PHYSICS->mcConfig().AddHadronicId(20443);
381 PHYSICS->mcConfig().AddHadronicId(20513);
382 PHYSICS->mcConfig().AddHadronicId(20523);
383 PHYSICS->mcConfig().AddHadronicId(20533);
384 PHYSICS->mcConfig().AddHadronicId(20543);
385 PHYSICS->mcConfig().AddHadronicId(20553);
386 PHYSICS->mcConfig().AddHadronicId(100443);
387 PHYSICS->mcConfig().AddHadronicId(100553);
388 PHYSICS->mcConfig().AddHadronicId(9900440);
389 PHYSICS->mcConfig().AddHadronicId(9900441);
390 PHYSICS->mcConfig().AddHadronicId(9900443);
391 PHYSICS->mcConfig().AddHadronicId(9900551);
392 PHYSICS->mcConfig().AddHadronicId(9900553);
393 PHYSICS->mcConfig().AddHadronicId(9910441);
394 PHYSICS->mcConfig().AddHadronicId(9910551);
395
396 // definition of the multiparticle "invisible"
397 PHYSICS->mcConfig().AddInvisibleId(-16);
398 PHYSICS->mcConfig().AddInvisibleId(-14);
399 PHYSICS->mcConfig().AddInvisibleId(-12);
400 PHYSICS->mcConfig().AddInvisibleId(12);
401 PHYSICS->mcConfig().AddInvisibleId(14);
402 PHYSICS->mcConfig().AddInvisibleId(16);
403 PHYSICS->mcConfig().AddInvisibleId(1000022);
404 PHYSICS->mcConfig().AddInvisibleId(1000039);
405
406 // Initializing PhysicsService for RECO
407 PHYSICS->recConfig().Reset();
408
409 // ========================= //
410 // ===== Signal region ===== //
411 // ========================= //
412
413 Manager()->AddRegionSelection("SRA");
414 Manager()->AddRegionSelection("SRA_L");
415 Manager()->AddRegionSelection("SRA_M");
416 Manager()->AddRegionSelection("SRA_H");
417
418 Manager()->AddRegionSelection("SRB");
419
420 Manager()->AddRegionSelection("SRC");
421 Manager()->AddRegionSelection("SRC_22");
422 Manager()->AddRegionSelection("SRC_24");
423 Manager()->AddRegionSelection("SRC_26");
424 Manager()->AddRegionSelection("SRC_28");
425
426 std::string SRA[] = {"SRA","SRA_L","SRA_M","SRA_H"};
427 std::string SRC[] = {"SRC","SRC_22","SRC_24","SRC_26","SRC_28"};
428 std::string SRAB[] = {"SRA","SRA_L","SRA_M","SRA_H","SRB"};
429 std::string SRAB1[] = {"SRA","SRB"};
430 std::string SRAC[] = {"SRA","SRA_L","SRA_M","SRA_H",
431 "SRC","SRC_22","SRC_24","SRC_26","SRC_28"};
432
433 // ====================== //
434 // ===== Selections ===== //
435 // ====================== //
436
437 // Common selections
438 Manager()->AddCut("$N_{lep} = 0$");
439
440
441 Manager()->AddCut("$N_{j} \\geq 6$", SRA);
442 Manager()->AddCut("$N_{j} \\geq 5$", "SRB");
443 Manager()->AddCut("$N_{j} \\geq 4$", SRC);
444
445 Manager()->AddCut("$N_{b} \\geq 4$", SRAB);
446 Manager()->AddCut("$N_{b} \\geq 3$", SRC);
447
448 Manager()->AddCut("$\\slashed{E}_T > 350$ [GeV]", SRAB);
449 Manager()->AddCut("$\\slashed{E}_T > 250$ [GeV]", SRC);
450
451 Manager()->AddCut("$min(\\Delta\\phi(j_{1-4},\\slashed{E}_T))>0.4$ [rad]");
452
453 // SRC selections
454 Manager()->AddCut("$\\slashed{E}_T$ Sig. $>22\\ [\\sqrt{\\rm GeV}]$", "SRC");
455 Manager()->AddCut("$\\slashed{E}_T$ Sig. $\\in [22,24]\\ [\\sqrt{\\rm GeV}]$", "SRC_22");
456 Manager()->AddCut("$\\slashed{E}_T$ Sig. $\\in [24,26]\\ [\\sqrt{\\rm GeV}]$", "SRC_24");
457 Manager()->AddCut("$\\slashed{E}_T$ Sig. $\\in [26,28]\\ [\\sqrt{\\rm GeV}]$", "SRC_26");
458 Manager()->AddCut("$\\slashed{E}_T$ Sig. $>28\\ [\\sqrt{\\rm GeV}]$", "SRC_28");
459
460 Manager()->AddCut("$\\tau^h$ veto", SRAB);
461
462 // SRA selections
463 Manager()->AddCut("$p^{b_1}_T > 200$ [GeV]", SRA);
464 Manager()->AddCut("$\\Delta R_{max}(b,b)>2.5$", SRA);
465 Manager()->AddCut("$\\Delta R_{max-min}(b,b)<2.5$", SRA);
466 Manager()->AddCut("$m(h_{cand})>80$ [GeV]", SRA);
467
468 // SRB selections
469 Manager()->AddCut("$m(h_{cand1},h_{cand2})_{avg}\\in[75,175]$ [GeV]", "SRB");
470 Manager()->AddCut("Leading jet non-b-tagged", "SRB");
471 Manager()->AddCut("$p^{j_1}_T > 350$ [GeV]", "SRB");
472 Manager()->AddCut("$|\\Delta\\phi(j_{1},\\slashed{E}_T)|>2.8$ [rad]", "SRB");
473
474 Manager()->AddCut("$m_{eff} > 1$ [TeV]", SRAB1);
475
476 Manager()->AddCut("$m_{eff} \\in [1,1.5]$ [TeV]", "SRA_L");
477 Manager()->AddCut("$m_{eff} \\in [1.5,2]$ [TeV]", "SRA_M");
478 Manager()->AddCut("$m_{eff} > 2$ [TeV]", "SRA_H");
479
480
481 // ====================== //
482 // ===== Histograms ===== //
483 // ====================== //
484
485 // Digitized histograms can be found in HEPData
486 // https://www.hepdata.net/record/ins1748602
487
488 Manager()->AddHisto("SRA_Meff", 11,800.0,3000., "SRA");
489 Manager()->AddHisto("SRA_Mh", 12,0.0,480., "SRA");
490
491 Manager()->AddHisto("SRB_PTj1", 9,50.0,950., "SRB");
492 Manager()->AddHisto("SRB_MhAvg", 16,50.0,450., "SRB");
493
494 Manager()->AddHisto("SRC_MET", 13,200.0,1500., "SRC");
495 Manager()->AddHisto("SRC_Sig", 19,17.0,36., "SRC");
496
497 return true;
498}
499
500// -----------------------------------------------------------------------------
501// Finalize
502// function called one time at the end of the analysis
503// -----------------------------------------------------------------------------
504void atlas_susy_2018_31::Finalize(const SampleFormat& summary,
505 const std::vector<SampleFormat>& files){}
506
507// -----------------------------------------------------------------------------
508// Execute
509// function called each time one event is read
510// -----------------------------------------------------------------------------
511bool atlas_susy_2018_31::Execute(SampleFormat& sample, const EventFormat& event)
512{
513 // Event weight
514 MAdouble64 EvWeight;
515 if(Configuration().IsNoEventWeight()) EvWeight=1.;
516 else if(event.mc()->weight()!=0.) EvWeight=event.mc()->weight();
517 else { return false;}
518 Manager()->InitializeForNewEvent(EvWeight);
519 // Empty event
520 if (event.rec()==0) {return true;}
521
522 std::vector<const RecJetFormat*> BaseJets, SignalBJets, SignalJets;
523 std::vector<const RecLeptonFormat*> preBaseMuons, preBaseElectrons;
524 std::vector<const RecLeptonFormat*> BaseMuons, BaseElectrons;
525 std::vector<const RecTauFormat*> SignalTaus;
526// std::vector<MALorentzVector> BaseObjects;
527
528 DEBUG << "============== New Event " << endmsg;
529 // Jets
530 for(MAuint32 ij=0; ij<event.rec()->jets().size(); ij++)
531 {
532 const RecJetFormat *CurrentJet = &(event.rec()->jets()[ij]);
533 if ( CurrentJet->pt() > 20.0 && CurrentJet->abseta() < 4.8)
534 {
535 BaseJets.push_back(CurrentJet);
536// if (CurrentJet->abseta() < 4.5)
537// BaseObjects.push_back(CurrentJet->momentum());
538 }
539 }
540 DEBUG << " * Reconstructed Jets..." << endmsg;
541
542 // Taus
543 // Candidates (tau) are identified as jets which have |eta| < 2.5 and
544 // less than five inner detector tracks of pT > 500 MeV.
545 for(MAuint32 ij=0; ij<event.rec()->taus().size(); ij++)
546 {
547 const RecTauFormat *CurrentTau = &(event.rec()->taus()[ij]);
548 if ( CurrentTau->pt() > 2.5 && CurrentTau->abseta()<2.5)
549 {
550 SignalTaus.push_back(CurrentTau);
551// if (CurrentTau->pt() > 10.)
552// BaseObjects.push_back(CurrentTau->momentum());
553 }
554 }
555 DEBUG << " * Reconstructed Taus..." << endmsg;
556
557 // Electron with pT > 4.5 GeV & eta < 2.47
558 for(MAuint32 ie=0; ie<event.rec()->electrons().size(); ie++)
559 {
560 const RecLeptonFormat *CurrentElectron = &(event.rec()->electrons()[ie]);
561 if( CurrentElectron->pt()>4.5 && CurrentElectron->abseta() < 2.47)
562 {
563 preBaseElectrons.push_back(CurrentElectron);
564// if (CurrentElectron->abseta() > 10.)
565// BaseObjects.push_back(CurrentElectron->momentum());
566 }
567 }
568 DEBUG << " * Reconstructed Electron..." << endmsg;
569
570 // Muons with pT > 4 GeV & eta < 2.5
571 for(MAuint32 im=0; im<event.rec()->muons().size(); im++)
572 {
573 const RecLeptonFormat *CurrentMuon = &(event.rec()->muons()[im]);
574 if( CurrentMuon->pt()>4. && CurrentMuon->abseta() < 2.5)
575 {
576 preBaseMuons.push_back(CurrentMuon);
577// if (CurrentMuon->abseta() > 10.)
578// BaseObjects.push_back(CurrentMuon->momentum());
579 }
580 }
581 DEBUG << " * Reconstructed Muon..." << endmsg;
582
583 // Lepton-jet overlap removal
584 BaseJets = Removal(BaseJets, preBaseElectrons, 0.2);
585 BaseElectrons = Removal(preBaseElectrons, BaseJets, 0.4);
586 BaseMuons = Removal(preBaseMuons, BaseJets, 0.4);
587
588 DEBUG << " * Removal..." << endmsg;
589
590 // Select signal jets and bjets from basejets
591 MAdouble64 HT = 0.0;
592 for (MAuint32 i=0; i<BaseJets.size(); i++)
593 {
594 if ( BaseJets[i]->pt() > 30.0 && BaseJets[i]->abseta() < 2.8)
595 {
596 SignalJets.push_back(BaseJets[i]);
597 HT += BaseJets[i]->pt();
598 if ( BaseJets[i]->abseta() < 2.5 && BaseJets[i]->btag())
599 SignalBJets.push_back(BaseJets[i]);
600 }
601 }
602 DEBUG << " * Based jets..." << endmsg;
603 DEBUG << " * Event reconstructed properly..." << endmsg;
604
605 // The inclusive MET & Meff
606 MALorentzVector pTmiss = event.rec()->MET().momentum();
607 MAdouble64 MET = pTmiss.Pt();
608 MAdouble64 Meff = HT + MET;
609
610 // Ordering everything
611 SORTER->sort(SignalJets, PTordering);
612 SORTER->sort(SignalBJets, PTordering);
613
614 MAint32 nBaselineLeptons = BaseMuons.size() + BaseElectrons.size();
615 MAuint32 nJets = SignalJets.size();
616 MAint32 nbJets = SignalBJets.size();
617
618 //==================//
619 //==== Cut Flow ====//
620 //==================//
621
622 DEBUG << " * Starting cut-flow..." << endmsg;
623
624 //if (MET<=35.) Manager()->SetCurrentEventWeight(0.);
625 if (MET<=250.) Manager()->SetCurrentEventWeight(EvWeight*0.80);
626 // ATL-DAQ-PUB-2016-001 p13 -> 35 GeV
627// if (!Manager()->ApplyCut(MET > 35.0, "$\\slashed{E}_T$ trigger")) return true;
628 if (!Manager()->ApplyCut(nBaselineLeptons == 0, "$N_{lep} = 0$")) return true;
629
630 if (!Manager()->ApplyCut(nJets >= 6, "$N_{j} \\geq 6$")) return true;
631 if (!Manager()->ApplyCut(nJets >= 5, "$N_{j} \\geq 5$")) return true;
632 if (!Manager()->ApplyCut(nJets >= 4, "$N_{j} \\geq 4$")) return true;
633
634 if (!Manager()->ApplyCut(nbJets >= 4, "$N_{b} \\geq 4$")) return true;
635 if (!Manager()->ApplyCut(nbJets >= 3, "$N_{b} \\geq 3$")) return true;
636
637 if (!Manager()->ApplyCut(MET > 350.0, "$\\slashed{E}_T > 350$ [GeV]")) return true;
638 if (!Manager()->ApplyCut(MET > 250.0, "$\\slashed{E}_T > 250$ [GeV]")) return true;
639
640 //minDPhis
641 MAdouble64 minDPhiJMET_4 = 99.;
642 MAdouble64 DPhiJMET_1 = 99.;
643 for (MAuint32 i=0; i<4; i++)
644 {
645 if (nJets >= i+1)
646 {
647 MAdouble64 DPhi = SignalJets[i]->dphi_0_pi(pTmiss);
648 if (DPhi < minDPhiJMET_4) minDPhiJMET_4 = DPhi;
649 }
650 }
651 if (nJets > 0) DPhiJMET_1 = SignalJets[0]->dphi_0_pi(pTmiss);
652
653 MAdouble64 minDPhiTauMet = 99.;
654 for (MAuint32 i=0; i<SignalTaus.size(); i++)
655 {
656 MAdouble64 DPhi = SignalTaus[i]->dphi_0_pi(pTmiss);
657 if (DPhi < minDPhiTauMet) minDPhiTauMet = DPhi;
658 }
659
660 if (!Manager()->ApplyCut(minDPhiJMET_4 > 0.4,
661 "$min(\\Delta\\phi(j_{1-4},\\slashed{E}_T))>0.4$ [rad]")) return true;
662
663
664 // ================================ //
665 // ======= MET Significance ======= //
666 // ================================ //
667
668 // see ref. https://cds.cern.ch/record/2630948
669
670// MAMatrix cov_sum(2,2), particle_u(2,2), particle_u_rot(2,2);
671// MALorentzVector softVec = pTmiss;
672// for (MAuint32 i=0; i<BaseObjects.size(); i++)
673// {
674// softVec += BaseObjects[i];
675// // needs object resolution for phi and pt
676// MAdouble64 pt_res=1., phi_res=1.;
677// particle_u[0][0] = std::pow(pt_res*BaseObjects[i].Pt(), 2);
678// particle_u[1][1] = std::pow(phi_res*BaseObjects[i].Phi(),2);
679// rotateXY(particle_u,particle_u_rot,BaseObjects[i].DeltaPhi(pTmiss));
680// add2x2Matrices(cov_sum,particle_u_rot);
681// }
682//
683// //add soft term resolution (fixed 10 GeV)
684// particle_u[0][0]=10*10;
685// particle_u[1][1]=10*10;
686// rotateXY(particle_u,particle_u_rot,pTmiss.DeltaPhi(softVec));
687// add2x2Matrices(cov_sum,particle_u_rot);
688//
689// //calculate significance
690// MAdouble64 varL = cov_sum[0][0];
691// MAdouble64 varT = cov_sum[1][1];
692// MAdouble64 covLT = cov_sum[0][1];
693// MAdouble64 significance = 0;
694// MAdouble64 rho = 0;
695// if( varL != 0 )
696// {
697// rho = covLT / sqrt( varL * varT ) ;
698// if (fabs( rho ) >= 0.9 ) rho=0; //too large - ignore it
699// significance = MET/sqrt((varL*(1-pow(rho,2))));
700// }
701//
702// INFO << "==============" << "\n"
703// << "varL = " << varL << "\n"
704// << "sig = " << significance << "\n"
705// << "rho = " << rho << endmsg;
706
707 MAdouble64 METtoHT = -1.;
708 if (nJets > 0.) METtoHT = MET/sqrt(HT);
709 Manager()->FillHisto("SRC_Sig", METtoHT);
710
711 /*
712 Manager()->FillHisto("SRC_Sig", significance);
713 if (!Manager()->ApplyCut(significance> 22., "$\\slashed{E}_T$ Sig. $>22\\ [\\sqrt{\\rm GeV}]$")) return true;
714 if (!Manager()->ApplyCut(significance>=22. && significance<=24., "$\\slashed{E}_T$ Sig. $\\in [22,24]\\ [\\sqrt{\\rm GeV}]$")) return true;
715 if (!Manager()->ApplyCut(significance>=24. && significance<=26., "$\\slashed{E}_T$ Sig. $\\in [24,26]\\ [\\sqrt{\\rm GeV}]$")) return true;
716 if (!Manager()->ApplyCut(significance>=26. && significance<=28., "$\\slashed{E}_T$ Sig. $\\in [26,28]\\ [\\sqrt{\\rm GeV}]$")) return true;
717 if (!Manager()->ApplyCut(significance> 28., "$\\slashed{E}_T$ Sig. $>28\\ [\\sqrt{\\rm GeV}]$")) return true;
718 */
719
720 //Tuned values for MET Sig approximation defined as METtoHT (Taken from ATLAS' public code in HEPData):
721 // See https://www.hepdata.net/record/ins1748602
722
723 if (!Manager()->ApplyCut(METtoHT> 18.4, "$\\slashed{E}_T$ Sig. $>22\\ [\\sqrt{\\rm GeV}]$")) return true;
724 if (!Manager()->ApplyCut(METtoHT>=18.4 && METtoHT<=21.0, "$\\slashed{E}_T$ Sig. $\\in [22,24]\\ [\\sqrt{\\rm GeV}]$")) return true;
725 if (!Manager()->ApplyCut(METtoHT>=21.0 && METtoHT<=23.55, "$\\slashed{E}_T$ Sig. $\\in [24,26]\\ [\\sqrt{\\rm GeV}]$")) return true;
726 if (!Manager()->ApplyCut(METtoHT>=23.55 && METtoHT<=26.1, "$\\slashed{E}_T$ Sig. $\\in [26,28]\\ [\\sqrt{\\rm GeV}]$")) return true;
727 if (!Manager()->ApplyCut(METtoHT> 27., "$\\slashed{E}_T$ Sig. $>28\\ [\\sqrt{\\rm GeV}]$")) return true;
728
729 Manager()->FillHisto("SRC_MET", MET);
730
731 if (!Manager()->ApplyCut(minDPhiTauMet>1.0471976, "$\\tau^h$ veto")) return true;
732 if (!Manager()->ApplyCut(SignalBJets[0]->pt() > 200.0, "$p^{b_1}_T > 200$ [GeV]")) return true;
733
734 // =============================== //
735 // ===== Max - Min Algorithm ===== //
736 // =============================== //
737
738 MAdouble64 maxDR = 0.;
739 MAdouble64 maxmin_M = -99.;
740 MAdouble64 max_M = -99.;
741 MAdouble64 max_avg_M = -99.;
742 MAuint32 imax = -1;
743 MAuint32 jmax = -1;
744 MAdouble64 maxminDR = 99.;
745 MAdouble64 maxmaxDR = -1.;
746 MAuint32 imin = -1;
747 MAuint32 jmin = -1;
748 MAuint32 imaxmax = -1;
749 MAuint32 jmaxmax = -1;
750
751 DEBUG << " * Starting Min-Max Algorithm..." << endmsg;
752
753 for(MAuint32 i=0; i<SignalBJets.size(); i++)
754 {
755 std::vector<MAdouble64> temp;
756 for(MAuint32 j=0; j<=i; j++)
757 {
758 MAdouble64 DR = SignalBJets[i]->dr(SignalBJets[j]);
759 if (i != j)
760 {
761 if (DR>maxDR)
762 {
763 maxDR = DR;
764 imax = i;
765 jmax = j;
766 }
767 } else {continue;}
768 }
769 }
770
771 max_M = (SignalBJets[imax]->momentum()+SignalBJets[jmax]->momentum()).M();
772
773 for(MAuint32 i=0; i<SignalBJets.size(); i++)
774 {
775 for(MAuint32 j=0; j<=i; j++)
776 {
777 if (i != j)
778 {
779 MAdouble64 DR = SignalBJets[i]->dr(SignalBJets[j]);
780 if (i != imax && j != jmax && DR > maxmaxDR)
781 {
782 imaxmax = i;
783 jmaxmax = j;
784 }
785 if (i != imax && j != jmax && DR < maxminDR)
786 {
787 maxminDR = DR;
788 imin = i;
789 jmin = j;
790 }
791 } else {continue;}
792 }
793 }
794
795 MAdouble64 maxmax_M = 0.;
796
797 if (SignalBJets.size() >= imaxmax && SignalBJets.size() >= jmaxmax && imaxmax >= 0 && jmaxmax >= 0)
798 maxmax_M = (SignalBJets[imaxmax]->momentum() + SignalBJets[jmaxmax]->momentum()).M();
799
800 if (max_M > 0. && maxmax_M > 0.)
801 max_avg_M = (max_M+maxmax_M)/2.;
802
803 if (SignalBJets.size() >= imin && SignalBJets.size() >= jmin && imin >= 0 && jmin >= 0)
804 maxmin_M = (SignalBJets[imin]->momentum() + SignalBJets[jmin]->momentum()).M();
805
806 if (!Manager()->ApplyCut(maxDR > 2.5, "$\\Delta R_{max}(b,b)>2.5$")) return true;
807 if (!Manager()->ApplyCut(maxminDR < 2.5, "$\\Delta R_{max-min}(b,b)<2.5$")) return true;
808
809 Manager()->FillHisto("SRA_Mh", maxmin_M);
810 if (!Manager()->ApplyCut(maxmin_M > 80.0, "$m(h_{cand})>80$ [GeV]")) return true;
811
812 Manager()->FillHisto("SRB_MhAvg", max_avg_M);
813 if (!Manager()->ApplyCut(max_avg_M <= 175.0 && max_avg_M >= 75.0,
814 "$m(h_{cand1},h_{cand2})_{avg}\\in[75,175]$ [GeV]")) return true;
815 if (!Manager()->ApplyCut(!SignalJets[0]->btag(), "Leading jet non-b-tagged")) return true;
816
817 Manager()->FillHisto("SRB_PTj1", SignalJets[0]->pt());
818 if (!Manager()->ApplyCut(SignalJets[0]->pt()>350.0, "$p^{j_1}_T > 350$ [GeV]")) return true;
819 if (!Manager()->ApplyCut(DPhiJMET_1>2.8, "$|\\Delta\\phi(j_{1},\\slashed{E}_T)|>2.8$ [rad]")) return true;
820
821 Manager()->FillHisto("SRA_Meff", Meff);
822 if (!Manager()->ApplyCut(Meff> 1000.0, "$m_{eff} > 1$ [TeV]")) return true;
823 if (!Manager()->ApplyCut(Meff>=1000.0 && Meff<=1500.0, "$m_{eff} \\in [1,1.5]$ [TeV]")) return true;
824 if (!Manager()->ApplyCut(Meff>=1500.0 && Meff<=2000.0, "$m_{eff} \\in [1.5,2]$ [TeV]")) return true;
825 if (!Manager()->ApplyCut(Meff> 2000.0, "$m_{eff} > 2$ [TeV]")) return true;
826
827 DEBUG << " * All cuts applied properly..." << endmsg;
828
829 return true;
830}
831
832// Overlap Removal
833template<typename T1, typename T2> std::vector<const T1*>
834 Removal(std::vector<const T1*> &v1, std::vector<const T2*> &v2,
835 const MAdouble64 &drmin)
836{
837 // Determining with objects should be removed
838 std::vector<bool> mask(v1.size(),false);
839 for (MAuint32 j=0;j<v1.size();j++)
840 for (MAuint32 i=0;i<v2.size();i++)
841 if (v2[i]->dr(v1[j]) < drmin)
842 {
843 mask[j]=true;
844 break;
845 }
846
847 // Building the cleaned container
848 std::vector<const T1*> cleaned_v1;
849 for (MAuint32 i=0;i<v1.size();i++)
850 if (!mask[i]) cleaned_v1.push_back(v1[i]);
851
852 return cleaned_v1;
853}
854
855// For MET significance calculation
856// static void rotateXY(MAMatrix &mat, MAMatrix &mat_new, MAdouble64 phi)
857// {
858// MAdouble64 c = cos(phi);
859// MAdouble64 s = sin(phi);
860// MAdouble64 cc = c*c;
861// MAdouble64 ss = s*s;
862// MAdouble64 cs = c*s;
863//
864// mat_new[0][0] = mat[0][0]*cc + mat[1][1]*ss - cs*(mat[1][0] + mat[0][1]);
865// mat_new[0][1] = mat[0][1]*cc - mat[1][0]*ss + cs*(mat[0][0] - mat[1][1]);
866// mat_new[1][0] = mat[1][0]*cc - mat[0][1]*ss + cs*(mat[0][0] - mat[1][1]);
867// mat_new[1][1] = mat[0][0]*ss + mat[1][1]*cc + cs*(mat[1][0] + mat[0][1]);
868// }
869//
870// static void add2x2Matrices(MAMatrix &mat1, MAMatrix &mat2)
871// {
872// for (MAuint32 ix=0; ix<2; ix++)
873// for (MAuint32 jx=0; jx<2; jx++)
874// mat1[ix][jx] += mat2[ix][jx];
875// }