rivet is hosted by Hepforge, IPPP Durham
Rivet 4.0.0
ParticleBaseUtils.hh
1#ifndef RIVET_PARTICLEBASEUTILS_HH
2#define RIVET_PARTICLEBASEUTILS_HH
3
4#include "Rivet/Tools/Utils.hh"
5#include "Rivet/ParticleBase.hh"
6
7namespace Rivet {
8
9
12
18
20 using ParticleBaseSelector = function<bool(const ParticleBase&)>;
22 using ParticleBaseSorter = function<bool(const ParticleBase&, const ParticleBase&)>;
23
24
27 virtual bool operator()(const ParticleBase& p) const = 0;
28 virtual ~BoolParticleBaseFunctor() {}
29 };
30
31
34 PtGtr(double pt) : ptcut(pt) { }
35 PtGtr(const FourMomentum& p) : ptcut(p.pT()) { }
36 bool operator()(const ParticleBase& p) const { return p.pT() > ptcut; }
37 double ptcut;
38 };
39 using pTGtr = PtGtr;
40 using ptGtr = PtGtr;
41
44 PtLess(const FourMomentum& p) : ptcut(p.pT()) { }
45 PtLess(double pt) : ptcut(pt) { }
46 bool operator()(const ParticleBase& p) const { return p.pT() < ptcut; }
47 double ptcut;
48 };
49 using pTLess = PtLess;
50 using ptLess = PtLess;
51
54 PtInRange(pair<double, double> ptcuts) : ptcut(ptcuts) { }
55 PtInRange(double ptlow, double pthigh) : PtInRange(make_pair(ptlow, pthigh)) { }
56 PtInRange(const FourMomentum& p1, const FourMomentum& p2) : PtInRange(p1.pT(), p2.pT()) { }
57 bool operator()(const ParticleBase& p) const { return p.pT() >= ptcut.first && p.pT() < ptcut.second; }
58 pair<double,double> ptcut;
59 };
60 using pTInRange = PtInRange;
61 using ptInRange = PtInRange;
62
63
66 EtaGtr(double eta) : etacut(eta) { }
67 EtaGtr(const FourMomentum& p) : etacut(p.eta()) { }
68 bool operator()(const ParticleBase& p) const { return p.eta() > etacut; }
69 double etacut;
70 };
71 using etaGtr = EtaGtr;
72
75 EtaLess(double eta) : etacut(eta) { }
76 EtaLess(const FourMomentum& p) : etacut(p.eta()) { }
77 bool operator()(const ParticleBase& p) const { return p.eta() < etacut; }
78 double etacut;
79 };
80 using etaLess = EtaLess;
81
84 EtaInRange(pair<double, double> etacuts) : etacut(etacuts) { }
85 EtaInRange(double etalow, double etahigh) : EtaInRange(make_pair(etalow, etahigh)) { }
86 EtaInRange(const FourMomentum& p1, const FourMomentum& p2) : EtaInRange(p1.eta(), p2.eta()) { }
87 bool operator()(const ParticleBase& p) const { return p.eta() >= etacut.first && p.eta() < etacut.second; }
88 pair<double,double> etacut;
89 };
90 using etaInRange = EtaInRange;
91
92
95 AbsEtaGtr(double abseta) : absetacut(abseta) { }
96 AbsEtaGtr(const FourMomentum& p) : absetacut(p.abseta()) { }
97 bool operator()(const ParticleBase& p) const { return p.abseta() > absetacut; }
98 double absetacut;
99 };
100 using absEtaGtr = AbsEtaGtr;
101 using absetaGtr = AbsEtaGtr;
102
105 AbsEtaLess(double abseta) : absetacut(abseta) { }
106 AbsEtaLess(const FourMomentum& p) : absetacut(p.abseta()) { }
107 bool operator()(const ParticleBase& p) const { return p.abseta() < absetacut; }
108 double absetacut;
109 };
110 using absEtaLess = AbsEtaLess;
111 using absetaLess = AbsEtaLess;
112
115 AbsEtaInRange(const pair<double, double>& absetacuts) : absetacut(absetacuts) { }
116 AbsEtaInRange(double absetalow, double absetahigh) : AbsEtaInRange(make_pair(absetalow, absetahigh)) { }
117 AbsEtaInRange(const FourMomentum& p1, const FourMomentum& p2) : AbsEtaInRange(p1.abseta(), p2.abseta()) { }
118 bool operator()(const ParticleBase& p) const { return p.abseta() >= absetacut.first && p.abseta() < absetacut.second; }
119 pair<double,double> absetacut;
120 };
123
124
127 RapGtr(double rap) : rapcut(rap) { }
128 RapGtr(const FourMomentum& p) : rapcut(p.rap()) { }
129 bool operator()(const ParticleBase& p) const { return p.rap() > rapcut; }
130 double rapcut;
131 };
132 using rapGtr = RapGtr;
133
136 RapLess(double rap) : rapcut(rap) { }
137 RapLess(const FourMomentum& p) : rapcut(p.rap()) { }
138 bool operator()(const ParticleBase& p) const { return p.rap() < rapcut; }
139 double rapcut;
140 };
141 using rapLess = RapLess;
142
145 RapInRange(const pair<double, double>& rapcuts) : rapcut(rapcuts) { }
146 RapInRange(double raplow, double raphigh) : RapInRange(make_pair(raplow, raphigh)) { }
147 RapInRange(const FourMomentum& p1, const FourMomentum& p2) : RapInRange(p1.rap(), p2.rap()) { }
148 bool operator()(const ParticleBase& p) const { return p.rap() >= rapcut.first && p.rap() < rapcut.second; }
149 pair<double,double> rapcut;
150 };
151 using rapInRange = RapInRange;
152
153
156 AbsRapGtr(double absrap) : absrapcut(absrap) { }
157 AbsRapGtr(const FourMomentum& p) : absrapcut(p.absrap()) { }
158 bool operator()(const ParticleBase& p) const { return p.absrap() > absrapcut; }
159 double absrapcut;
160 };
161 using absRapGtr = AbsRapGtr;
162 using absrapGtr = AbsRapGtr;
163
166 AbsRapLess(double absrap) : absrapcut(absrap) { }
167 AbsRapLess(const FourMomentum& p) : absrapcut(p.absrap()) { }
168 bool operator()(const ParticleBase& p) const { return p.absrap() < absrapcut; }
169 double absrapcut;
170 };
171 using absRapLess = AbsRapLess;
172 using absrapLess = AbsRapLess;
173
176 AbsRapInRange(const pair<double, double>& absrapcuts) : absrapcut(absrapcuts) { }
177 AbsRapInRange(double absraplow, double absraphigh) : AbsRapInRange(make_pair(absraplow, absraphigh)) { }
178 AbsRapInRange(const FourMomentum& p1, const FourMomentum& p2) : AbsRapInRange(p1.absrap(), p2.absrap()) { }
179 bool operator()(const ParticleBase& p) const { return p.absrap() >= absrapcut.first && p.absrap() < absrapcut.second; }
180 pair<double,double> absrapcut;
181 };
184
185
187
188
191 DeltaRGtr(const ParticleBase& vec, double dr, RapScheme scheme=PSEUDORAPIDITY)
192 : refvec(vec.mom()), drcut(dr), rapscheme(scheme) { }
193 DeltaRGtr(const FourMomentum& vec, double dr, RapScheme scheme=PSEUDORAPIDITY)
194 : refvec(vec), drcut(dr), rapscheme(scheme) { }
195 DeltaRGtr(const Vector3& vec, double dr)
196 : drcut(dr), rapscheme(PSEUDORAPIDITY) { refvec.setPx(vec.x()); refvec.setPy(vec.y()); refvec.setPz(vec.z()); }
197 bool operator()(const ParticleBase& p) const { return deltaR(p, refvec, rapscheme) > drcut; }
198 FourMomentum refvec;
199 double drcut;
200 RapScheme rapscheme;
201 };
202 using deltaRGtr = DeltaRGtr;
203
206 DeltaRLess(const ParticleBase& vec, double dr, RapScheme scheme=PSEUDORAPIDITY)
207 : refvec(vec.mom()), drcut(dr), rapscheme(scheme) { }
208 DeltaRLess(const FourMomentum& vec, double dr, RapScheme scheme=PSEUDORAPIDITY)
209 : refvec(vec), drcut(dr), rapscheme(scheme) { }
210 DeltaRLess(const Vector3& vec, double dr)
211 : drcut(dr), rapscheme(PSEUDORAPIDITY) { refvec.setPx(vec.x()); refvec.setPy(vec.y()); refvec.setPz(vec.z()); }
212 bool operator()(const ParticleBase& p) const { return deltaR(p, refvec, rapscheme) < drcut; }
213 FourMomentum refvec;
214 double drcut;
215 RapScheme rapscheme;
216 };
217 using deltaRLess = DeltaRLess;
218
221 DeltaRInRange(const ParticleBase& vec, const pair<double,double>& dr, RapScheme scheme=PSEUDORAPIDITY)
222 : refvec(vec.mom()), drcut(dr), rapscheme(scheme) { }
223 DeltaRInRange(const ParticleBase& vec, double drmin, double drmax, RapScheme scheme=PSEUDORAPIDITY)
224 : DeltaRInRange(vec, make_pair(drmin, drmax), scheme) { }
225 DeltaRInRange(const FourMomentum& vec, const pair<double,double>& dr, RapScheme scheme=PSEUDORAPIDITY)
226 : refvec(vec), drcut(dr), rapscheme(scheme) { }
227 DeltaRInRange(const FourMomentum& vec, double drmin, double drmax, RapScheme scheme=PSEUDORAPIDITY)
228 : DeltaRInRange(vec, make_pair(drmin, drmax), scheme) { }
229 DeltaRInRange(const Vector3& vec, const pair<double,double>& dr)
230 : drcut(dr), rapscheme(PSEUDORAPIDITY) { refvec.setPx(vec.x()); refvec.setPy(vec.y()); refvec.setPz(vec.z()); }
231 DeltaRInRange(const Vector3& vec, double drmin, double drmax)
232 : DeltaRInRange(vec, make_pair(drmin, drmax)) { }
233 bool operator()(const ParticleBase& p) const {
234 const double dR = deltaR(p, refvec, rapscheme);
235 return dR >= drcut.first && dR < drcut.second;
236 }
237 FourMomentum refvec;
238 pair<double,double> drcut;
239 RapScheme rapscheme;
240 };
242
243
246 DeltaPhiGtr(const ParticleBase& vec, double dphi)
247 : refvec(vec.p3()), dphicut(dphi) { }
248 DeltaPhiGtr(const FourMomentum& vec, double dphi)
249 : refvec(vec.p3()), dphicut(dphi) { }
250 DeltaPhiGtr(const Vector3& vec, double dphi)
251 : refvec(vec), dphicut(dphi) { }
252 bool operator()(const ParticleBase& p) const { return deltaPhi(p, refvec) > dphicut; }
253 Vector3 refvec;
254 double dphicut;
255 };
256 using deltaPhiGtr = DeltaPhiGtr;
257
260 DeltaPhiLess(const ParticleBase& vec, double dphi)
261 : refvec(vec.p3()), dphicut(dphi) { }
262 DeltaPhiLess(const FourMomentum& vec, double dphi)
263 : refvec(vec.p3()), dphicut(dphi) { }
264 DeltaPhiLess(const Vector3& vec, double dphi)
265 : refvec(vec), dphicut(dphi) { }
266 bool operator()(const ParticleBase& p) const { return deltaPhi(p, refvec) < dphicut; }
267 Vector3 refvec;
268 double dphicut;
269 };
271
274 DeltaPhiInRange(const ParticleBase& vec, const pair<double,double>& dphi)
275 : refvec(vec.mom()), dphicut(dphi) { }
276 DeltaPhiInRange(const ParticleBase& vec, double dphimin, double dphimax)
277 : DeltaPhiInRange(vec, make_pair(dphimin, dphimax)) { }
278 DeltaPhiInRange(const FourMomentum& vec, const pair<double,double>& dphi)
279 : refvec(vec), dphicut(dphi) { }
280 DeltaPhiInRange(const FourMomentum& vec, double dphimin, double dphimax)
281 : DeltaPhiInRange(vec, make_pair(dphimin, dphimax)) { }
282 DeltaPhiInRange(const Vector3& vec, const pair<double,double>& dphi)
283 : refvec(vec), dphicut(dphi) { }
284 DeltaPhiInRange(const Vector3& vec, double dphimin, double dphimax)
285 : DeltaPhiInRange(vec, make_pair(dphimin, dphimax)) { }
286 bool operator()(const ParticleBase& p) const {
287 const double dphi = deltaPhi(p, refvec);
288 return dphi >= dphicut.first && dphi < dphicut.second;
289 }
290 Vector3 refvec;
291 pair<double,double> dphicut;
292 };
294
295
298 DeltaEtaGtr(const ParticleBase& vec, double deta)
299 : refvec(vec.p3()), detacut(deta) { }
300 DeltaEtaGtr(const FourMomentum& vec, double deta)
301 : refvec(vec.p3()), detacut(deta) { }
302 DeltaEtaGtr(const Vector3& vec, double deta)
303 : refvec(vec), detacut(deta) { }
304 bool operator()(const ParticleBase& p) const { return std::abs(deltaEta(p, refvec)) > detacut; }
305 Vector3 refvec;
306 double detacut;
307 };
308 using deltaEtaGtr = DeltaEtaGtr;
309
312 DeltaEtaLess(const ParticleBase& vec, double deta)
313 : refvec(vec.p3()), detacut(deta) { }
314 DeltaEtaLess(const FourMomentum& vec, double deta)
315 : refvec(vec.p3()), detacut(deta) { }
316 DeltaEtaLess(const Vector3& vec, double deta)
317 : refvec(vec), detacut(deta) { }
318 bool operator()(const ParticleBase& p) const { return std::abs(deltaEta(p, refvec)) < detacut; }
319 Vector3 refvec;
320 double detacut;
321 };
323
326 DeltaEtaInRange(const ParticleBase& vec, const pair<double,double>& deta)
327 : refvec(vec.mom()), detacut(deta) { }
328 DeltaEtaInRange(const ParticleBase& vec, double detamin, double detamax)
329 : DeltaEtaInRange(vec, make_pair(detamin, detamax)) { }
330 DeltaEtaInRange(const FourMomentum& vec, const pair<double,double>& deta)
331 : refvec(vec), detacut(deta) { }
332 DeltaEtaInRange(const FourMomentum& vec, double detamin, double detamax)
333 : DeltaEtaInRange(vec, make_pair(detamin, detamax)) { }
334 DeltaEtaInRange(const Vector3& vec, const pair<double,double>& deta)
335 : refvec(vec), detacut(deta) { }
336 DeltaEtaInRange(const Vector3& vec, double detamin, double detamax)
337 : DeltaEtaInRange(vec, make_pair(detamin, detamax)) { }
338 bool operator()(const ParticleBase& p) const {
339 const double deta = deltaEta(p, refvec);
340 return deta >= detacut.first && deta < detacut.second;
341 }
342 Vector3 refvec;
343 pair<double,double> detacut;
344 };
346
347
350 DeltaRapGtr(const ParticleBase& vec, double drap)
351 : refvec(vec.mom()), drapcut(drap) { }
352 DeltaRapGtr(const FourMomentum& vec, double drap)
353 : refvec(vec), drapcut(drap) { }
354 bool operator()(const ParticleBase& p) const { return std::abs(deltaRap(p, refvec)) > drapcut; }
355 FourMomentum refvec;
356 double drapcut;
357 };
358 using deltaRapGtr = DeltaRapGtr;
359
362 DeltaRapLess(const ParticleBase& vec, double drap)
363 : refvec(vec.mom()), drapcut(drap) { }
364 DeltaRapLess(const FourMomentum& vec, double drap)
365 : refvec(vec), drapcut(drap) { }
366 bool operator()(const ParticleBase& p) const { return std::abs(deltaRap(p, refvec)) < drapcut; }
367 FourMomentum refvec;
368 double drapcut;
369 };
371
374 DeltaRapInRange(const ParticleBase& vec, const pair<double,double>& drap)
375 : refvec(vec.mom()), drapcut(drap) { }
376 DeltaRapInRange(const ParticleBase& vec, double drapmin, double drapmax)
377 : DeltaRapInRange(vec, make_pair(drapmin, drapmax)) { }
378 DeltaRapInRange(const FourMomentum& vec, const pair<double,double>& drap)
379 : refvec(vec), drapcut(drap) { }
380 DeltaRapInRange(const FourMomentum& vec, double drapmin, double drapmax)
381 : DeltaRapInRange(vec, make_pair(drapmin, drapmax)) { }
382 bool operator()(const ParticleBase& p) const {
383 const double drap = deltaRap(p, refvec);
384 return drap >= drapcut.first && drap < drapcut.second;
385 }
386 FourMomentum refvec;
387 pair<double,double> drapcut;
388 };
390
392
393
399
402 virtual double operator()(const ParticleBase& p) const = 0;
403 virtual ~DoubleParticleBaseFunctor() {}
404 };
405
408 DeltaRWRT(const ParticleBase& pb, RapScheme scheme=PSEUDORAPIDITY) : p(pb.mom()), rapscheme(scheme) {}
409 DeltaRWRT(const FourMomentum& p4, RapScheme scheme=PSEUDORAPIDITY) : p(p4), rapscheme(scheme) {}
410 DeltaRWRT(const Vector3& p3) : p(p3.mod(), p3.x(), p3.y(), p3.z()), rapscheme(PSEUDORAPIDITY) {}
411 double operator()(const ParticleBase& pb) const { return deltaR(p, pb, rapscheme); }
412 double operator()(const FourMomentum& p4) const { return deltaR(p, p4, rapscheme); }
413 double operator()(const Vector3& p3) const { return deltaR(p, p3); }
414 const FourMomentum p;
415 RapScheme rapscheme;
416 };
417 using deltaRWRT = DeltaRWRT;
418
421 DeltaPhiWRT(const ParticleBase& pb) : p(pb.mom().vector3()) {}
422 DeltaPhiWRT(const FourMomentum& p4) : p(p4.vector3()) {}
423 DeltaPhiWRT(const Vector3& p3) : p(p3) {}
424 double operator()(const ParticleBase& pb) const { return deltaPhi(p, pb); }
425 double operator()(const FourMomentum& p4) const { return deltaPhi(p, p4); }
426 double operator()(const Vector3& p3) const { return deltaPhi(p, p3); }
427 const Vector3 p;
428 };
429 using deltaPhiWRT = DeltaPhiWRT;
430
433 DeltaEtaWRT(const ParticleBase& pb) : p(pb.mom().vector3()) {}
434 DeltaEtaWRT(const FourMomentum& p4) : p(p4.vector3()) {}
435 DeltaEtaWRT(const Vector3& p3) : p(p3) {}
436 double operator()(const ParticleBase& pb) const { return deltaEta(p, pb); }
437 double operator()(const FourMomentum& p4) const { return deltaEta(p, p4); }
438 double operator()(const Vector3& p3) const { return deltaEta(p, p3); }
439 const Vector3 p;
440 };
441 using deltaEtaWRT = DeltaEtaWRT;
442
445 AbsDeltaEtaWRT(const ParticleBase& pb) : p(pb.mom().vector3()) {}
446 AbsDeltaEtaWRT(const FourMomentum& p4) : p(p4.vector3()) {}
447 AbsDeltaEtaWRT(const Vector3& p3) : p(p3) {}
448 double operator()(const ParticleBase& pb) const { return fabs(deltaEta(p, pb)); }
449 double operator()(const FourMomentum& p4) const { return fabs(deltaEta(p, p4)); }
450 double operator()(const Vector3& p3) const { return fabs(deltaEta(p, p3)); }
451 const Vector3 p;
452 };
454
457 DeltaRapWRT(const ParticleBase& pb) : p(pb.mom()) {}
458 DeltaRapWRT(const FourMomentum& p4) : p(p4) {}
459 double operator()(const ParticleBase& pb) const { return deltaRap(p, pb); }
460 double operator()(const FourMomentum& p4) const { return deltaRap(p, p4); }
461 const FourMomentum p;
462 };
463 using deltaRapWRT = DeltaRapWRT;
464
467 AbsDeltaRapWRT(const ParticleBase& pb) : p(pb.mom()) {}
468 AbsDeltaRapWRT(const FourMomentum& p4) : p(p4) {}
469 double operator()(const ParticleBase& pb) const { return fabs(deltaRap(p, pb)); }
470 double operator()(const FourMomentum& p4) const { return fabs(deltaRap(p, p4)); }
471 const FourMomentum p;
472 };
474
476
477
480
481 template<typename PBCONTAINER1, typename PBCONTAINER2>
482 inline void idiscardIfAny(PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare,
483 typename std::function<bool(const typename PBCONTAINER1::value_type&,
484 const typename PBCONTAINER2::value_type&)> fn) {
485 for (const auto& pbcmp : tocompare) {
486 idiscard(tofilter, [&](const typename PBCONTAINER1::value_type& pbfilt){ return fn(pbfilt, pbcmp); });
487 }
488 }
489
490 template<typename PBCONTAINER1, typename PBCONTAINER2>
491 inline PBCONTAINER1 discardIfAny(const PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare,
492 typename std::function<bool(const typename PBCONTAINER1::value_type&,
493 const typename PBCONTAINER2::value_type&)> fn) {
494 PBCONTAINER1 tmp{tofilter};
495 idiscardIfAny(tmp, tocompare, fn);
496 return tmp;
497 }
498
499
500 template<typename PBCONTAINER1, typename PBCONTAINER2>
501 inline PBCONTAINER1 selectIfAny(const PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare,
502 typename std::function<bool(const typename PBCONTAINER1::value_type&,
503 const typename PBCONTAINER2::value_type&)> fn) {
504 PBCONTAINER1 selected;
505 for (const auto& pbfilt : tofilter) {
506 if (any(tocompare, [&](const typename PBCONTAINER2::value_type& pbcmp){ return fn(pbfilt, pbcmp); })) {
507 selected += pbfilt;
508 }
509 }
510 return selected;
511 }
512
513 template<typename PBCONTAINER1, typename PBCONTAINER2>
514 inline void iselectIfAny(PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare,
515 typename std::function<bool(const typename PBCONTAINER1::value_type&,
516 const typename PBCONTAINER2::value_type&)> fn) {
517 tofilter = selectIfAny(tofilter, tocompare, fn);
518 }
519
520
521
522 template<typename PBCONTAINER1, typename PBCONTAINER2>
523 inline PBCONTAINER1 discardIfAll(const PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare,
524 typename std::function<bool(const typename PBCONTAINER1::value_type&,
525 const typename PBCONTAINER2::value_type&)> fn) {
526 PBCONTAINER1 selected;
527 for (const auto& pbfilt : tofilter) {
528 if (!all(tocompare, [&](const typename PBCONTAINER2::value_type& pbcmp){ return fn(pbfilt, pbcmp); })) {
529 selected += pbfilt;
530 }
531 }
532 return selected;
533 }
534
535 template<typename PBCONTAINER1, typename PBCONTAINER2>
536 inline void idiscardIfAll(PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare,
537 typename std::function<bool(const typename PBCONTAINER1::value_type&,
538 const typename PBCONTAINER2::value_type&)> fn) {
539 tofilter = discardIfAll(tofilter, tocompare, fn);
540 }
541
542
543 template<typename PBCONTAINER1, typename PBCONTAINER2>
544 inline PBCONTAINER1 selectIfAll(const PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare,
545 typename std::function<bool(const typename PBCONTAINER1::value_type&,
546 const typename PBCONTAINER2::value_type&)> fn) {
547 PBCONTAINER1 selected;
548 for (const auto& pbfilt : tofilter) {
549 if (all(tocompare, [&](const typename PBCONTAINER2::value_type& pbcmp){ return fn(pbfilt, pbcmp); })) {
550 selected += pbfilt;
551 }
552 }
553 return selected;
554 }
555
556 template<typename PBCONTAINER1, typename PBCONTAINER2>
557 inline void iselectIfAll(PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare,
558 typename std::function<bool(const typename PBCONTAINER1::value_type&,
559 const typename PBCONTAINER2::value_type&)> fn) {
560 tofilter = selectIfAll(tofilter, tocompare, fn);
561 }
562
564
565
568
569 template<typename PBCONTAINER1, typename PBCONTAINER2>
570 inline void idiscardIfAnyDeltaRLess(PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare, double dR) {
571 for (const typename PBCONTAINER2::value_type& pb : tocompare) {
572 idiscard(tofilter, deltaRLess(pb, dR));
573 }
574 }
575
576 template<typename PBCONTAINER1, typename PBCONTAINER2>
577 inline PBCONTAINER1 discardIfAnyDeltaRLess(const PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare, double dR) {
578 PBCONTAINER1 tmp{tofilter};
579 idiscardIfAnyDeltaRLess(tmp, tocompare, dR);
580 return tmp;
581 }
582
583 template<typename PBCONTAINER1, typename PBCONTAINER2>
584 inline void idiscardIfAnyDeltaPhiLess(PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare, double dphi) {
585 for (const typename PBCONTAINER2::value_type& pb : tocompare) {
586 idiscard(tofilter, deltaPhiLess(pb, dphi));
587 }
588 }
589
590 template<typename PBCONTAINER1, typename PBCONTAINER2>
591 inline PBCONTAINER1 discardIfAnyDeltaPhiLess(const PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare, double dphi) {
592 PBCONTAINER1 tmp{tofilter};
593 idiscardIfAnyDeltaPhiLess(tmp, tocompare, dphi);
594 return tmp;
595 }
596
597
598
599 template<typename PBCONTAINER1, typename PBCONTAINER2>
600 inline PBCONTAINER1 selectIfAnyDeltaRLess(const PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare, double dR) {
601 PBCONTAINER1 selected;
602 for (const typename PBCONTAINER1::value_type& f : tofilter) {
603 if (any(tocompare, deltaRLess(f, dR))) selected.push_back(f);
604 }
605 return selected;
606 }
607
608 template<typename PBCONTAINER1, typename PBCONTAINER2>
609 inline void iselectIfAnyDeltaRLess(PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare, double dR) {
610 tofilter = selectIfAnyDeltaRLess(tofilter, tocompare, dR);
611 }
612
613
614 template<typename PBCONTAINER1, typename PBCONTAINER2>
615 inline PBCONTAINER1 selectIfAnyDeltaPhiLess(const PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare, double dphi) {
616 PBCONTAINER1 selected;
617 for (const typename PBCONTAINER1::value_type& f : tofilter) {
618 if (any(tocompare, deltaPhiLess(f, dphi))) selected.push_back(f);
619 }
620 return selected;
621 }
622
623 template<typename PBCONTAINER1, typename PBCONTAINER2>
624 inline void iselectIfAnyDeltaPhiLess(PBCONTAINER1& tofilter, const PBCONTAINER2& tocompare, double dphi) {
625 tofilter = selectIfAnyDeltaPhiLess(tofilter, tocompare, dphi);
626 }
627
628
630
632
633
639 namespace Kin {
640
642 inline FourMomentum mom(const ParticleBase& p) { return p.mom(); }
644 inline FourMomentum p4(const ParticleBase& p) { return p.mom(); }
645
647 inline Vector3 p3(const ParticleBase& p) { return p.p3(); }
648
650 inline Vector3 pTvec(const ParticleBase& p) { return p.pTvec(); }
651
653 inline double p(const ParticleBase& p) { return p.p(); }
654
656 inline double pT(const ParticleBase& p) { return p.pT(); }
657
659 inline double Et(const ParticleBase& p) { return p.Et(); }
660
662 inline double eta(const ParticleBase& p) { return p.eta(); }
663
665 inline double abseta(const ParticleBase& p) { return p.abseta(); }
666
668 inline double rap(const ParticleBase& p) { return p.rap(); }
669
671 inline double absrap(const ParticleBase& p) { return p.absrap(); }
672
674 inline double mass(const ParticleBase& p) { return p.mass(); }
675
676
677 // /// Unbound function access to pair pT
678 // inline double pairPt(const ParticleBase& p1, const ParticleBase& p2) { return (p1.mom() + p2.mom()).pT(); }
679
680 // /// Unbound function access to pair mass
681 // inline double pairMass(const ParticleBase& p1, const ParticleBase& p2) { return (p1.mom() + p2.mom()).mass(); }
682
683
685 inline double mass(const ParticleBase& p, const FourMomentum& p4) {
686 return mass(p.mom(), p4);
687 }
688
690 inline double mass(const FourMomentum& p4, const ParticleBase& p) {
691 return mass(p4, p.mom());
692 }
693
695 inline double mass(const ParticleBase& p1, const ParticleBase& p2) {
696 return mass(p1.mom(), p2.mom());
697 }
698
700 inline double mass2(const ParticleBase& p, const P4& p4) {
701 return mass2(p.mom(), p4);
702 }
703
705 inline double mass2(const P4& p4, const ParticleBase& p) {
706 return mass2(p4, p.mom());
707 }
708
710 inline double mass2(const ParticleBase& p1, const ParticleBase& p2) {
711 return mass2(p1.mom(), p2.mom());
712 }
713
718 inline double mT(const ParticleBase& p, const P4& p4) {
719 return mT(p.mom(), p4);
720 }
721
726 inline double mT(const P4& p4, const ParticleBase& p) {
727 return mT(p4, p.mom());
728 }
729
734 inline double mT(const ParticleBase& p1, const ParticleBase& p2) {
735 return mT(p1.mom(), p2.mom());
736 }
737
739 inline double pT(const ParticleBase& p, const P4& p4) {
740 return pT(p.mom(), p4);
741 }
742
744 inline double pT(const P4& p4, const ParticleBase& p) {
745 return pT(p4, p.mom());
746 }
747
749 inline double pT(const ParticleBase& p1, const ParticleBase& p2) {
750 return pT(p1.mom(), p2.mom());
751 }
752
753 }
754
755 // Import Kin namespace into Rivet
756 using namespace Kin;
757
759
760
763
768 // auto closestMassIndex = std::bind(closestMatchIndex, std::placeholders::_1, Kin::mass, std::placeholders::_2, std::placeholders::_3);
769 template <typename CONTAINER, typename = isCIterable<CONTAINER>>
770 inline int closestMassIndex(CONTAINER&& c, double mtarget, double mmin=-DBL_MAX, double mmax=DBL_MAX) {
771 return closestMatchIndex(std::forward<CONTAINER>(c), Kin::mass, mtarget, mmin, mmax);
772 }
773
778 template <typename CONTAINER1, typename CONTAINER2, typename = isCIterable<CONTAINER1, CONTAINER2>>
779 inline pair<int,int> closestMassIndices(CONTAINER1&& c1, CONTAINER2&& c2,
780 double mtarget, double mmin=-DBL_MAX, double mmax=DBL_MAX) {
781 return closestMatchIndices(std::forward<CONTAINER1>(c1),
782 std::forward<CONTAINER2>(c2), Kin::mass, mtarget, mmin, mmax);
783 }
784
789 template <typename CONTAINER, typename T, typename = isCIterable<CONTAINER>>
790 inline int closestMassIndex(CONTAINER&& c, const T& x,
791 double mtarget, double mmin=-DBL_MAX, double mmax=DBL_MAX) {
792 using FN = double(const ParticleBase&,const T&);
793 return closestMatchIndex<CONTAINER, T, FN>(std::forward<CONTAINER>(c), x, Kin::mass, mtarget, mmin, mmax);
794 }
795
796
801 template <typename CONTAINER, typename T, typename = isCIterable<CONTAINER>>
802 inline int closestMassIndex(T&& x, CONTAINER&& c,
803 double mtarget, double mmin=-DBL_MAX, double mmax=DBL_MAX) {
804 return closestMatchIndex(std::forward<T>(x), std::forward<CONTAINER>(c), Kin::mass, mtarget, mmin, mmax);
805 }
806
808
809
811
812}
813
814#endif
Specialized version of the FourVector with momentum/energy functionality.
Definition Vector4.hh:316
FourMomentum & setPz(double pz)
Set z-component of momentum .
Definition Vector4.hh:370
FourMomentum & setPy(double py)
Set y-component of momentum .
Definition Vector4.hh:364
double rap() const
Alias for rapidity.
Definition Vector4.hh:611
FourMomentum & setPx(double px)
Set x-component of momentum .
Definition Vector4.hh:358
double absrap() const
Absolute rapidity.
Definition Vector4.hh:620
double pT() const
Calculate the transverse momentum .
Definition Vector4.hh:643
Vector3 vector3() const
Get the spatial part of the 4-vector as a 3-vector.
Definition Vector4.hh:177
double abseta() const
Get the directly (alias).
Definition Vector4.hh:174
double eta() const
Synonym for pseudorapidity.
Definition Vector4.hh:167
Base class for particle-like things like Particle and Jet.
Definition ParticleBase.hh:13
const FourMomentum & mom() const
Get equivalent single momentum four-vector (const) (alias).
Definition ParticleBase.hh:39
ThreeMomentum p3() const
Get the 3-momentum directly.
Definition ParticleBase.hh:108
Three-dimensional specialisation of Vector.
Definition Vector3.hh:40
double mod() const
Calculate the modulus of a vector. .
Definition VectorN.hh:95
pair< int, int > closestMatchIndices(CONTAINER1 &&c1, CONTAINER2 &&c2, FN &&fn, double target, double minval=-DBL_MAX, double maxval=DBL_MAX)
Return the indices from two vectors which best match fn(c1[i], c2[j]) to the target value.
Definition Utils.hh:788
int closestMatchIndex(CONTAINER &&c, FN &&fn, double target, double minval=-DBL_MAX, double maxval=DBL_MAX)
Return the index from a vector which best matches fn(c[i]) to the target value.
Definition Utils.hh:759
bool any(const CONTAINER &c)
Return true if x is true for any x in container c, otherwise false.
Definition Utils.hh:330
bool all(const CONTAINER &c)
Return true if x is true for all x in container c, otherwise false.
Definition Utils.hh:355
Jets & idiscard(Jets &jets, const Cut &c)
Filter a jet collection in-place to the subset that fails the supplied Cut.
double mass2(const FourMomentum &a, const FourMomentum &b)
Calculate mass^2 of two 4-vectors.
Definition Vector4.hh:1466
function< bool(const ParticleBase &)> ParticleBaseSelector
std::function instantiation for functors taking a ParticleBase and returning a bool
Definition ParticleBaseUtils.hh:20
function< bool(const ParticleBase &, const ParticleBase &)> ParticleBaseSorter
std::function instantiation for functors taking two ParticleBase and returning a bool
Definition ParticleBaseUtils.hh:22
pair< int, int > closestMassIndices(CONTAINER1 &&c1, CONTAINER2 &&c2, double mtarget, double mmin=-DBL_MAX, double mmax=DBL_MAX)
Return the indices from two vectors which best match fn(c1[i], c2[j]) to the target value.
Definition ParticleBaseUtils.hh:779
int closestMassIndex(CONTAINER &&c, double mtarget, double mmin=-DBL_MAX, double mmax=DBL_MAX)
Return the index from a vector which best matches mass(c[i]) to the target value.
Definition ParticleBaseUtils.hh:770
Definition MC_CENT_PPB_Projections.hh:10
double deltaR(double rap1, double phi1, double rap2, double phi2)
Definition MathUtils.hh:698
double deltaPhi(double phi1, double phi2, bool sign=false)
Calculate the difference between two angles in radians.
Definition MathUtils.hh:668
double deltaEta(double eta1, double eta2, bool sign=false)
Definition MathUtils.hh:676
double mT(double pT1, double pT2, double dphi)
Definition MathUtils.hh:720
RapScheme
Enum for rapidity variable to be used in calculating , applying rapidity cuts, etc.
Definition MathConstants.hh:46
double deltaRap(double y1, double y2, bool sign=false)
Definition MathUtils.hh:684
Calculator of with respect to a given momentum.
Definition ParticleBaseUtils.hh:444
Calculator of with respect to a given momentum.
Definition ParticleBaseUtils.hh:466
Abs pseudorapidity greater-than functor.
Definition ParticleBaseUtils.hh:94
Abs pseudorapidity in-range functor.
Definition ParticleBaseUtils.hh:114
Abs pseudorapidity momentum less-than functor.
Definition ParticleBaseUtils.hh:104
Abs rapidity greater-than functor.
Definition ParticleBaseUtils.hh:155
Abs rapidity in-range functor.
Definition ParticleBaseUtils.hh:175
Abs rapidity momentum less-than functor.
Definition ParticleBaseUtils.hh:165
Base type for Particle -> bool functors.
Definition ParticleBaseUtils.hh:26
(with respect to another momentum, vec) greater-than functor
Definition ParticleBaseUtils.hh:297
(with respect to another 4-momentum, vec) in-range functor
Definition ParticleBaseUtils.hh:325
(with respect to another momentum, vec) less-than functor
Definition ParticleBaseUtils.hh:311
Calculator of with respect to a given momentum.
Definition ParticleBaseUtils.hh:432
(with respect to another momentum, vec) greater-than functor
Definition ParticleBaseUtils.hh:245
(with respect to another 4-momentum, vec) in-range functor
Definition ParticleBaseUtils.hh:273
(with respect to another momentum, vec) less-than functor
Definition ParticleBaseUtils.hh:259
Calculator of with respect to a given momentum.
Definition ParticleBaseUtils.hh:420
(with respect to another 4-momentum, vec) greater-than functor
Definition ParticleBaseUtils.hh:190
(with respect to another 4-momentum, vec) in-range functor
Definition ParticleBaseUtils.hh:220
(with respect to another 4-momentum, vec) less-than functor
Definition ParticleBaseUtils.hh:205
Calculator of with respect to a given momentum.
Definition ParticleBaseUtils.hh:407
(with respect to another momentum, vec) greater-than functor
Definition ParticleBaseUtils.hh:349
(with respect to another 4-momentum, vec) in-range functor
Definition ParticleBaseUtils.hh:373
(with respect to another momentum, vec) less-than functor
Definition ParticleBaseUtils.hh:361
Calculator of with respect to a given momentum.
Definition ParticleBaseUtils.hh:456
Base type for Particle -> double functors.
Definition ParticleBaseUtils.hh:401
Pseudorapidity greater-than functor.
Definition ParticleBaseUtils.hh:65
Pseudorapidity in-range functor.
Definition ParticleBaseUtils.hh:83
Pseudorapidity less-than functor.
Definition ParticleBaseUtils.hh:74
Transverse momentum greater-than functor.
Definition ParticleBaseUtils.hh:33
Transverse momentum in-range functor.
Definition ParticleBaseUtils.hh:53
Transverse momentum less-than functor.
Definition ParticleBaseUtils.hh:43
Rapidity greater-than functor.
Definition ParticleBaseUtils.hh:126
Rapidity in-range functor.
Definition ParticleBaseUtils.hh:144
Rapidity momentum less-than functor.
Definition ParticleBaseUtils.hh:135