Fork me on GitHub

source: svn/trunk/external/fastjet/Selector.cc@ 1352

Last change on this file since 1352 was 1332, checked in by Pavel Demin, 11 years ago

upgrade fastjet to 3.0.6

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id Revision Date
File size: 47.0 KB
Line 
1//STARTHEADER
2// $Id: Selector.cc 1332 2013-11-20 20:52:59Z pavel $
3//
4// Copyright (c) 2005-2011, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
5//
6//----------------------------------------------------------------------
7// This file is part of FastJet.
8//
9// FastJet is free software; you can redistribute it and/or modify
10// it under the terms of the GNU General Public License as published by
11// the Free Software Foundation; either version 2 of the License, or
12// (at your option) any later version.
13//
14// The algorithms that underlie FastJet have required considerable
15// development and are described in hep-ph/0512210. If you use
16// FastJet as part of work towards a scientific publication, please
17// include a citation to the FastJet paper.
18//
19// FastJet is distributed in the hope that it will be useful,
20// but WITHOUT ANY WARRANTY; without even the implied warranty of
21// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
22// GNU General Public License for more details.
23//
24// You should have received a copy of the GNU General Public License
25// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
26//----------------------------------------------------------------------
27//ENDHEADER
28
29
30#include <sstream>
31#include <algorithm>
32#include "fastjet/Selector.hh"
33#ifndef __FJCORE__
34#include "fastjet/GhostedAreaSpec.hh" // for area support
35#endif // __FJCORE__
36
37using namespace std;
38
39FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
40
41//----------------------------------------------------------------------
42// implementations of some of the more complex bits of Selector
43//----------------------------------------------------------------------
44
45// implementation of the operator() acting on a vector of jets
46std::vector<PseudoJet> Selector::operator()(const std::vector<PseudoJet> & jets) const {
47 std::vector<PseudoJet> result;
48 const SelectorWorker * worker_local = validated_worker();
49 if (worker_local->applies_jet_by_jet()) {
50 //if (false) {
51 // for workers that apply jet by jet, this is more efficient
52 for (std::vector<PseudoJet>::const_iterator jet = jets.begin();
53 jet != jets.end(); jet++) {
54 if (worker_local->pass(*jet)) result.push_back(*jet);
55 }
56 } else {
57 // for workers that can only be applied to entire vectors,
58 // go through the following
59 std::vector<const PseudoJet *> jetptrs(jets.size());
60 for (unsigned i = 0; i < jets.size(); i++) {
61 jetptrs[i] = & jets[i];
62 }
63 worker_local->terminator(jetptrs);
64 for (unsigned i = 0; i < jetptrs.size(); i++) {
65 if (jetptrs[i]) result.push_back(jets[i]);
66 }
67 }
68 return result;
69}
70
71
72//----------------------------------------------------------------------
73// count the number of jets that pass the cuts
74unsigned int Selector::count(const std::vector<PseudoJet> & jets) const {
75 unsigned n = 0;
76 const SelectorWorker * worker_local = validated_worker();
77
78 // separate strategies according to whether the worker applies jet by jet
79 if (worker_local->applies_jet_by_jet()) {
80 for (unsigned i = 0; i < jets.size(); i++) {
81 if (worker_local->pass(jets[i])) n++;
82 }
83 } else {
84 std::vector<const PseudoJet *> jetptrs(jets.size());
85 for (unsigned i = 0; i < jets.size(); i++) {
86 jetptrs[i] = & jets[i];
87 }
88 worker_local->terminator(jetptrs);
89 for (unsigned i = 0; i < jetptrs.size(); i++) {
90 if (jetptrs[i]) n++;
91 }
92 }
93
94 return n;
95}
96
97
98//----------------------------------------------------------------------
99// sift the input jets into two vectors -- those that pass the selector
100// and those that do not
101void Selector::sift(const std::vector<PseudoJet> & jets,
102 std::vector<PseudoJet> & jets_that_pass,
103 std::vector<PseudoJet> & jets_that_fail
104 ) const {
105 const SelectorWorker * worker_local = validated_worker();
106
107 jets_that_pass.clear();
108 jets_that_fail.clear();
109
110 // separate strategies according to whether the worker applies jet by jet
111 if (worker_local->applies_jet_by_jet()) {
112 for (unsigned i = 0; i < jets.size(); i++) {
113 if (worker_local->pass(jets[i])) {
114 jets_that_pass.push_back(jets[i]);
115 } else {
116 jets_that_fail.push_back(jets[i]);
117 }
118 }
119 } else {
120 std::vector<const PseudoJet *> jetptrs(jets.size());
121 for (unsigned i = 0; i < jets.size(); i++) {
122 jetptrs[i] = & jets[i];
123 }
124 worker_local->terminator(jetptrs);
125 for (unsigned i = 0; i < jetptrs.size(); i++) {
126 if (jetptrs[i]) {
127 jets_that_pass.push_back(jets[i]);
128 } else {
129 jets_that_fail.push_back(jets[i]);
130 }
131 }
132 }
133}
134
135#ifndef __FJCORE__
136// area using default ghost area
137double Selector::area() const{
138 return area(gas::def_ghost_area);
139}
140
141// implementation of the Selector's area function
142double Selector::area(double ghost_area) const{
143 if (! is_geometric()) throw InvalidArea();
144
145 // has area will already check we've got a valid worker
146 if (_worker->has_known_area()) return _worker->known_area();
147
148 // generate a set of "ghosts"
149 double rapmin, rapmax;
150 get_rapidity_extent(rapmin, rapmax);
151 GhostedAreaSpec ghost_spec(rapmin, rapmax, 1, ghost_area);
152 std::vector<PseudoJet> ghosts;
153 ghost_spec.add_ghosts(ghosts);
154
155 // check what passes the selection
156 return ghost_spec.ghost_area() * ((*this)(ghosts)).size();
157}
158#endif // __FJCORE__
159
160
161//----------------------------------------------------------------------
162// implementations of some of the more complex bits of SelectorWorker
163//----------------------------------------------------------------------
164// check if it has a finite area
165bool SelectorWorker::has_finite_area() const {
166 if (! is_geometric()) return false;
167 double rapmin, rapmax;
168 get_rapidity_extent(rapmin, rapmax);
169 return (rapmax != std::numeric_limits<double>::infinity())
170 && (-rapmin != std::numeric_limits<double>::infinity());
171}
172
173
174
175//----------------------------------------------------------------------
176// very basic set of selectors (at the moment just the identity!)
177//----------------------------------------------------------------------
178
179//----------------------------------------------------------------------
180/// helper for selecting the n hardest jets
181class SW_Identity : public SelectorWorker {
182public:
183 /// ctor with specification of the number of objects to keep
184 SW_Identity(){}
185
186 /// just let everything pass
187 virtual bool pass(const PseudoJet &) const {
188 return true;
189 }
190
191 /// For each jet that does not pass the cuts, this routine sets the
192 /// pointer to 0.
193 virtual void terminator(vector<const PseudoJet *> &) const {
194 // everything passes, hence nothing to nullify
195 return;
196 }
197
198 /// returns a description of the worker
199 virtual string description() const { return "Identity";}
200
201 /// strictly speaking, this is geometric
202 virtual bool is_geometric() const { return true;}
203};
204
205
206// returns an "identity" selector that lets everything pass
207Selector SelectorIdentity() {
208 return Selector(new SW_Identity);
209}
210
211
212//----------------------------------------------------------------------
213// selector and workers for operators
214//----------------------------------------------------------------------
215
216//----------------------------------------------------------------------
217/// helper for combining selectors with a logical not
218class SW_Not : public SelectorWorker {
219public:
220 /// ctor
221 SW_Not(const Selector & s) : _s(s) {}
222
223 /// return a copy of the current object
224 virtual SelectorWorker* copy(){ return new SW_Not(*this);}
225
226 /// returns true if a given object passes the selection criterium
227 /// this has to be overloaded by derived workers
228 virtual bool pass(const PseudoJet & jet) const {
229 // make sure that the "pass" can be applied on both selectors
230 if (!applies_jet_by_jet())
231 throw Error("Cannot apply this selector worker to an individual jet");
232
233 return ! _s.pass(jet);
234 }
235
236 /// returns true if this can be applied jet by jet
237 virtual bool applies_jet_by_jet() const {return _s.applies_jet_by_jet();}
238
239 /// select the jets in the list that pass both selectors
240 virtual void terminator(vector<const PseudoJet *> & jets) const {
241 // if we can apply the selector jet-by-jet, call the base selector
242 // that does exactly that
243 if (applies_jet_by_jet()){
244 SelectorWorker::terminator(jets);
245 return;
246 }
247
248 // check the effect of the selector we want to negate
249 vector<const PseudoJet *> s_jets = jets;
250 _s.worker()->terminator(s_jets);
251
252 // now apply the negation: all the jets that pass the base
253 // selector (i.e. are not NULL) have to be set to NULL
254 for (unsigned int i=0; i<s_jets.size(); i++){
255 if (s_jets[i]) jets[i] = NULL;
256 }
257 }
258
259 /// returns a description of the worker
260 virtual string description() const {
261 ostringstream ostr;
262 ostr << "!(" << _s.description() << ")";
263 return ostr.str();
264 }
265
266 /// is geometric if the underlying selector is
267 virtual bool is_geometric() const { return _s.is_geometric();}
268
269 /// returns true if the worker can be set_referenced
270 virtual bool takes_reference() const { return _s.takes_reference();}
271
272 /// set the reference jet for this selector
273 virtual void set_reference(const PseudoJet &ref) { _s.set_reference(ref);}
274
275protected:
276 Selector _s;
277};
278
279
280// logical not applied on a selector
281Selector operator!(const Selector & s) {
282 return Selector(new SW_Not(s));
283}
284
285
286//----------------------------------------------------------------------
287/// Base class for binary operators
288class SW_BinaryOperator: public SelectorWorker {
289public:
290 /// ctor
291 SW_BinaryOperator(const Selector & s1, const Selector & s2) : _s1(s1), _s2(s2) {
292 // stores info for more efficient access to the selector's properties
293
294 // we can apply jet by jet only if this is the case for both sub-selectors
295 _applies_jet_by_jet = _s1.applies_jet_by_jet() && _s2.applies_jet_by_jet();
296
297 // the selector takes a reference if either of the sub-selectors does
298 _takes_reference = _s1.takes_reference() || _s2.takes_reference();
299
300 // we have a well-defined area provided the two objects have one
301 _is_geometric = _s1.is_geometric() && _s2.is_geometric();
302 }
303
304 /// returns true if this can be applied jet by jet
305 virtual bool applies_jet_by_jet() const {return _applies_jet_by_jet;}
306
307 /// returns true if this takes a reference jet
308 virtual bool takes_reference() const{
309 return _takes_reference;
310 }
311
312 /// sets the reference jet
313 virtual void set_reference(const PseudoJet &centre){
314 _s1.set_reference(centre);
315 _s2.set_reference(centre);
316 }
317
318 /// check if it has a finite area
319 virtual bool is_geometric() const { return _is_geometric;}
320
321protected:
322 Selector _s1, _s2;
323 bool _applies_jet_by_jet;
324 bool _takes_reference;
325 bool _is_geometric;
326};
327
328
329
330//----------------------------------------------------------------------
331/// helper for combining selectors with a logical and
332class SW_And: public SW_BinaryOperator {
333public:
334 /// ctor
335 SW_And(const Selector & s1, const Selector & s2) : SW_BinaryOperator(s1,s2){}
336
337 /// return a copy of this
338 virtual SelectorWorker* copy(){ return new SW_And(*this);}
339
340 /// returns true if a given object passes the selection criterium
341 /// this has to be overloaded by derived workers
342 virtual bool pass(const PseudoJet & jet) const {
343 // make sure that the "pass" can be applied on both selectors
344 if (!applies_jet_by_jet())
345 throw Error("Cannot apply this selector worker to an individual jet");
346
347 return _s1.pass(jet) && _s2.pass(jet);
348 }
349
350 /// select the jets in the list that pass both selectors
351 virtual void terminator(vector<const PseudoJet *> & jets) const {
352 // if we can apply the selector jet-by-jet, call the base selector
353 // that does exactly that
354 if (applies_jet_by_jet()){
355 SelectorWorker::terminator(jets);
356 return;
357 }
358
359 // check the effect of the first selector
360 vector<const PseudoJet *> s1_jets = jets;
361 _s1.worker()->terminator(s1_jets);
362
363 // apply the second
364 _s2.worker()->terminator(jets);
365
366 // terminate the jets that wiuld be terminated by _s1
367 for (unsigned int i=0; i<jets.size(); i++){
368 if (! s1_jets[i]) jets[i] = NULL;
369 }
370 }
371
372 /// returns the rapidity range for which it may return "true"
373 virtual void get_rapidity_extent(double & rapmin, double & rapmax) const {
374 double s1min, s1max, s2min, s2max;
375 _s1.get_rapidity_extent(s1min, s1max);
376 _s2.get_rapidity_extent(s2min, s2max);
377 rapmax = min(s1max, s2max);
378 rapmin = max(s1min, s2min);
379 }
380
381 /// returns a description of the worker
382 virtual string description() const {
383 ostringstream ostr;
384 ostr << "(" << _s1.description() << " && " << _s2.description() << ")";
385 return ostr.str();
386 }
387};
388
389
390// logical and between two selectors
391Selector operator&&(const Selector & s1, const Selector & s2) {
392 return Selector(new SW_And(s1,s2));
393}
394
395
396
397//----------------------------------------------------------------------
398/// helper for combining selectors with a logical or
399class SW_Or: public SW_BinaryOperator {
400public:
401 /// ctor
402 SW_Or(const Selector & s1, const Selector & s2) : SW_BinaryOperator(s1,s2) {}
403
404 /// return a copy of this
405 virtual SelectorWorker* copy(){ return new SW_Or(*this);}
406
407 /// returns true if a given object passes the selection criterium
408 /// this has to be overloaded by derived workers
409 virtual bool pass(const PseudoJet & jet) const {
410 // make sure that the "pass" can be applied on both selectors
411 if (!applies_jet_by_jet())
412 throw Error("Cannot apply this selector worker to an individual jet");
413
414 return _s1.pass(jet) || _s2.pass(jet);
415 }
416
417 /// returns true if this can be applied jet by jet
418 virtual bool applies_jet_by_jet() const {
419 // watch out, even though it's the "OR" selector, to be applied jet
420 // by jet, both the base selectors need to be jet-by-jet-applicable,
421 // so the use of a && in the line below
422 return _s1.applies_jet_by_jet() && _s2.applies_jet_by_jet();
423 }
424
425 /// select the jets in the list that pass both selectors
426 virtual void terminator(vector<const PseudoJet *> & jets) const {
427 // if we can apply the selector jet-by-jet, call the base selector
428 // that does exactly that
429 if (applies_jet_by_jet()){
430 SelectorWorker::terminator(jets);
431 return;
432 }
433
434 // check the effect of the first selector
435 vector<const PseudoJet *> s1_jets = jets;
436 _s1.worker()->terminator(s1_jets);
437
438 // apply the second
439 _s2.worker()->terminator(jets);
440
441 // resurrect any jet that has been terminated by the second one
442 // and not by the first one
443 for (unsigned int i=0; i<jets.size(); i++){
444 if (s1_jets[i]) jets[i] = s1_jets[i];
445 }
446 }
447
448 /// returns a description of the worker
449 virtual string description() const {
450 ostringstream ostr;
451 ostr << "(" << _s1.description() << " || " << _s2.description() << ")";
452 return ostr.str();
453 }
454
455 /// returns the rapidity range for which it may return "true"
456 virtual void get_rapidity_extent(double & rapmin, double & rapmax) const {
457 double s1min, s1max, s2min, s2max;
458 _s1.get_rapidity_extent(s1min, s1max);
459 _s2.get_rapidity_extent(s2min, s2max);
460 rapmax = max(s1max, s2max);
461 rapmin = min(s1min, s2min);
462 }
463};
464
465
466// logical or between two selectors
467Selector operator ||(const Selector & s1, const Selector & s2) {
468 return Selector(new SW_Or(s1,s2));
469}
470
471//----------------------------------------------------------------------
472/// helper for multiplying two selectors (in an operator-like way)
473class SW_Mult: public SW_And {
474public:
475 /// ctor
476 SW_Mult(const Selector & s1, const Selector & s2) : SW_And(s1,s2) {}
477
478 /// return a copy of this
479 virtual SelectorWorker* copy(){ return new SW_Mult(*this);}
480
481 /// select the jets in the list that pass both selectors
482 virtual void terminator(vector<const PseudoJet *> & jets) const {
483 // if we can apply the selector jet-by-jet, call the base selector
484 // that does exactly that
485 if (applies_jet_by_jet()){
486 SelectorWorker::terminator(jets);
487 return;
488 }
489
490 // first apply _s2
491 _s2.worker()->terminator(jets);
492
493 // then apply _s1
494 _s1.worker()->terminator(jets);
495 }
496
497 /// returns a description of the worker
498 virtual string description() const {
499 ostringstream ostr;
500 ostr << "(" << _s1.description() << " * " << _s2.description() << ")";
501 return ostr.str();
502 }
503};
504
505
506// logical and between two selectors
507Selector operator*(const Selector & s1, const Selector & s2) {
508 return Selector(new SW_Mult(s1,s2));
509}
510
511
512//----------------------------------------------------------------------
513// selector and workers for kinematic cuts
514//----------------------------------------------------------------------
515
516//----------------------------------------------------------------------
517// a series of basic classes that allow easy implementations of
518// min, max and ranges on a quantity-to-be-defined
519
520// generic holder for a quantity
521class QuantityBase{
522public:
523 QuantityBase(double q) : _q(q){}
524 virtual ~QuantityBase(){}
525 virtual double operator()(const PseudoJet & jet ) const =0;
526 virtual string description() const =0;
527 virtual bool is_geometric() const { return false;}
528 virtual double comparison_value() const {return _q;}
529 virtual double description_value() const {return comparison_value();}
530protected:
531 double _q;
532};
533
534// generic holder for a squared quantity
535class QuantitySquareBase : public QuantityBase{
536public:
537 QuantitySquareBase(double sqrtq) : QuantityBase(sqrtq*sqrtq), _sqrtq(sqrtq){}
538 virtual double description_value() const {return _sqrtq;}
539protected:
540 double _sqrtq;
541};
542
543// generic_quantity >= minimum
544template<typename QuantityType>
545class SW_QuantityMin : public SelectorWorker{
546public:
547 /// detfault ctor (initialises the pt cut)
548 SW_QuantityMin(double qmin) : _qmin(qmin) {}
549
550 /// returns true is the given object passes the selection pt cut
551 virtual bool pass(const PseudoJet & jet) const {return _qmin(jet) >= _qmin.comparison_value();}
552
553 /// returns a description of the worker
554 virtual string description() const {
555 ostringstream ostr;
556 ostr << _qmin.description() << " >= " << _qmin.description_value();
557 return ostr.str();
558 }
559
560 virtual bool is_geometric() const { return _qmin.is_geometric();}
561
562protected:
563 QuantityType _qmin; ///< the cut
564};
565
566
567// generic_quantity <= maximum
568template<typename QuantityType>
569class SW_QuantityMax : public SelectorWorker {
570public:
571 /// detfault ctor (initialises the pt cut)
572 SW_QuantityMax(double qmax) : _qmax(qmax) {}
573
574 /// returns true is the given object passes the selection pt cut
575 virtual bool pass(const PseudoJet & jet) const {return _qmax(jet) <= _qmax.comparison_value();}
576
577 /// returns a description of the worker
578 virtual string description() const {
579 ostringstream ostr;
580 ostr << _qmax.description() << " <= " << _qmax.description_value();
581 return ostr.str();
582 }
583
584 virtual bool is_geometric() const { return _qmax.is_geometric();}
585
586protected:
587 QuantityType _qmax; ///< the cut
588};
589
590
591// generic quantity in [minimum:maximum]
592template<typename QuantityType>
593class SW_QuantityRange : public SelectorWorker {
594public:
595 /// detfault ctor (initialises the pt cut)
596 SW_QuantityRange(double qmin, double qmax) : _qmin(qmin), _qmax(qmax) {}
597
598 /// returns true is the given object passes the selection pt cut
599 virtual bool pass(const PseudoJet & jet) const {
600 double q = _qmin(jet); // we could identically use _qmax
601 return (q >= _qmin.comparison_value()) && (q <= _qmax.comparison_value());
602 }
603
604 /// returns a description of the worker
605 virtual string description() const {
606 ostringstream ostr;
607 ostr << _qmin.description_value() << " <= " << _qmin.description() << " <= " << _qmax.description_value();
608 return ostr.str();
609 }
610
611 virtual bool is_geometric() const { return _qmin.is_geometric();}
612
613protected:
614 QuantityType _qmin; // the lower cut
615 QuantityType _qmax; // the upper cut
616};
617
618
619//----------------------------------------------------------------------
620/// helper class for selecting on pt
621class QuantityPt2 : public QuantitySquareBase{
622public:
623 QuantityPt2(double pt) : QuantitySquareBase(pt){}
624 virtual double operator()(const PseudoJet & jet ) const { return jet.perp2();}
625 virtual string description() const {return "pt";}
626};
627
628// returns a selector for a minimum pt
629Selector SelectorPtMin(double ptmin) {
630 return Selector(new SW_QuantityMin<QuantityPt2>(ptmin));
631}
632
633// returns a selector for a maximum pt
634Selector SelectorPtMax(double ptmax) {
635 return Selector(new SW_QuantityMax<QuantityPt2>(ptmax));
636}
637
638// returns a selector for a pt range
639Selector SelectorPtRange(double ptmin, double ptmax) {
640 return Selector(new SW_QuantityRange<QuantityPt2>(ptmin, ptmax));
641}
642
643
644//----------------------------------------------------------------------
645/// helper class for selecting on transverse energy
646class QuantityEt2 : public QuantitySquareBase{
647public:
648 QuantityEt2(double Et) : QuantitySquareBase(Et){}
649 virtual double operator()(const PseudoJet & jet ) const { return jet.Et2();}
650 virtual string description() const {return "Et";}
651};
652
653// returns a selector for a minimum Et
654Selector SelectorEtMin(double Etmin) {
655 return Selector(new SW_QuantityMin<QuantityEt2>(Etmin));
656}
657
658// returns a selector for a maximum Et
659Selector SelectorEtMax(double Etmax) {
660 return Selector(new SW_QuantityMax<QuantityEt2>(Etmax));
661}
662
663// returns a selector for a Et range
664Selector SelectorEtRange(double Etmin, double Etmax) {
665 return Selector(new SW_QuantityRange<QuantityEt2>(Etmin, Etmax));
666}
667
668
669//----------------------------------------------------------------------
670/// helper class for selecting on energy
671class QuantityE : public QuantityBase{
672public:
673 QuantityE(double E) : QuantityBase(E){}
674 virtual double operator()(const PseudoJet & jet ) const { return jet.E();}
675 virtual string description() const {return "E";}
676};
677
678// returns a selector for a minimum E
679Selector SelectorEMin(double Emin) {
680 return Selector(new SW_QuantityMin<QuantityE>(Emin));
681}
682
683// returns a selector for a maximum E
684Selector SelectorEMax(double Emax) {
685 return Selector(new SW_QuantityMax<QuantityE>(Emax));
686}
687
688// returns a selector for a E range
689Selector SelectorERange(double Emin, double Emax) {
690 return Selector(new SW_QuantityRange<QuantityE>(Emin, Emax));
691}
692
693
694//----------------------------------------------------------------------
695/// helper class for selecting on mass
696class QuantityM2 : public QuantitySquareBase{
697public:
698 QuantityM2(double m) : QuantitySquareBase(m){}
699 virtual double operator()(const PseudoJet & jet ) const { return jet.m2();}
700 virtual string description() const {return "mass";}
701};
702
703// returns a selector for a minimum mass
704Selector SelectorMassMin(double mmin) {
705 return Selector(new SW_QuantityMin<QuantityM2>(mmin));
706}
707
708// returns a selector for a maximum mass
709Selector SelectorMassMax(double mmax) {
710 return Selector(new SW_QuantityMax<QuantityM2>(mmax));
711}
712
713// returns a selector for a mass range
714Selector SelectorMassRange(double mmin, double mmax) {
715 return Selector(new SW_QuantityRange<QuantityM2>(mmin, mmax));
716}
717
718
719
720//----------------------------------------------------------------------
721/// helper for selecting on rapidities: quantity
722class QuantityRap : public QuantityBase{
723public:
724 QuantityRap(double rap) : QuantityBase(rap){}
725 virtual double operator()(const PseudoJet & jet ) const { return jet.rap();}
726 virtual string description() const {return "rap";}
727 virtual bool is_geometric() const { return true;}
728};
729
730
731/// helper for selecting on rapidities: min
732class SW_RapMin : public SW_QuantityMin<QuantityRap>{
733public:
734 SW_RapMin(double rapmin) : SW_QuantityMin<QuantityRap>(rapmin){}
735 virtual void get_rapidity_extent(double &rapmin, double & rapmax) const{
736 rapmax = std::numeric_limits<double>::max();
737 rapmin = _qmin.comparison_value();
738 }
739};
740
741/// helper for selecting on rapidities: max
742class SW_RapMax : public SW_QuantityMax<QuantityRap>{
743public:
744 SW_RapMax(double rapmax) : SW_QuantityMax<QuantityRap>(rapmax){}
745 virtual void get_rapidity_extent(double &rapmin, double & rapmax) const{
746 rapmax = _qmax.comparison_value();
747 rapmin = -std::numeric_limits<double>::max();
748 }
749};
750
751/// helper for selecting on rapidities: range
752class SW_RapRange : public SW_QuantityRange<QuantityRap>{
753public:
754 SW_RapRange(double rapmin, double rapmax) : SW_QuantityRange<QuantityRap>(rapmin, rapmax){
755 assert(rapmin<=rapmax);
756 }
757 virtual void get_rapidity_extent(double &rapmin, double & rapmax) const{
758 rapmax = _qmax.comparison_value();
759 rapmin = _qmin.comparison_value();
760 }
761 virtual bool has_known_area() const { return true;} ///< the area is analytically known
762 virtual double known_area() const {
763 return twopi * (_qmax.comparison_value()-_qmin.comparison_value());
764 }
765};
766
767// returns a selector for a minimum rapidity
768Selector SelectorRapMin(double rapmin) {
769 return Selector(new SW_RapMin(rapmin));
770}
771
772// returns a selector for a maximum rapidity
773Selector SelectorRapMax(double rapmax) {
774 return Selector(new SW_RapMax(rapmax));
775}
776
777// returns a selector for a rapidity range
778Selector SelectorRapRange(double rapmin, double rapmax) {
779 return Selector(new SW_RapRange(rapmin, rapmax));
780}
781
782
783//----------------------------------------------------------------------
784/// helper for selecting on |rapidities|
785class QuantityAbsRap : public QuantityBase{
786public:
787 QuantityAbsRap(double absrap) : QuantityBase(absrap){}
788 virtual double operator()(const PseudoJet & jet ) const { return abs(jet.rap());}
789 virtual string description() const {return "|rap|";}
790 virtual bool is_geometric() const { return true;}
791};
792
793
794/// helper for selecting on |rapidities|: max
795class SW_AbsRapMax : public SW_QuantityMax<QuantityAbsRap>{
796public:
797 SW_AbsRapMax(double absrapmax) : SW_QuantityMax<QuantityAbsRap>(absrapmax){}
798 virtual void get_rapidity_extent(double &rapmin, double & rapmax) const{
799 rapmax = _qmax.comparison_value();
800 rapmin = -_qmax.comparison_value();
801 }
802 virtual bool has_known_area() const { return true;} ///< the area is analytically known
803 virtual double known_area() const {
804 return twopi * 2 * _qmax.comparison_value();
805 }
806};
807
808/// helper for selecting on |rapidities|: max
809class SW_AbsRapRange : public SW_QuantityRange<QuantityAbsRap>{
810public:
811 SW_AbsRapRange(double absrapmin, double absrapmax) : SW_QuantityRange<QuantityAbsRap>(absrapmin, absrapmax){}
812 virtual void get_rapidity_extent(double &rapmin, double & rapmax) const{
813 rapmax = _qmax.comparison_value();
814 rapmin = -_qmax.comparison_value();
815 }
816 virtual bool has_known_area() const { return true;} ///< the area is analytically known
817 virtual double known_area() const {
818 return twopi * 2 * (_qmax.comparison_value()-max(_qmin.comparison_value(),0.0)); // this should handle properly absrapmin<0
819 }
820};
821
822// returns a selector for a minimum |rapidity|
823Selector SelectorAbsRapMin(double absrapmin) {
824 return Selector(new SW_QuantityMin<QuantityAbsRap>(absrapmin));
825}
826
827// returns a selector for a maximum |rapidity|
828Selector SelectorAbsRapMax(double absrapmax) {
829 return Selector(new SW_AbsRapMax(absrapmax));
830}
831
832// returns a selector for a |rapidity| range
833Selector SelectorAbsRapRange(double rapmin, double rapmax) {
834 return Selector(new SW_AbsRapRange(rapmin, rapmax));
835}
836
837
838//----------------------------------------------------------------------
839/// helper for selecting on pseudo-rapidities
840class QuantityEta : public QuantityBase{
841public:
842 QuantityEta(double eta) : QuantityBase(eta){}
843 virtual double operator()(const PseudoJet & jet ) const { return jet.eta();}
844 virtual string description() const {return "eta";}
845 // virtual bool is_geometric() const { return true;} // not strictly only y and phi-dependent
846};
847
848// returns a selector for a pseudo-minimum rapidity
849Selector SelectorEtaMin(double etamin) {
850 return Selector(new SW_QuantityMin<QuantityEta>(etamin));
851}
852
853// returns a selector for a pseudo-maximum rapidity
854Selector SelectorEtaMax(double etamax) {
855 return Selector(new SW_QuantityMax<QuantityEta>(etamax));
856}
857
858// returns a selector for a pseudo-rapidity range
859Selector SelectorEtaRange(double etamin, double etamax) {
860 return Selector(new SW_QuantityRange<QuantityEta>(etamin, etamax));
861}
862
863
864//----------------------------------------------------------------------
865/// helper for selecting on |pseudo-rapidities|
866class QuantityAbsEta : public QuantityBase{
867public:
868 QuantityAbsEta(double abseta) : QuantityBase(abseta){}
869 virtual double operator()(const PseudoJet & jet ) const { return abs(jet.eta());}
870 virtual string description() const {return "|eta|";}
871 virtual bool is_geometric() const { return true;}
872};
873
874// returns a selector for a minimum |pseudo-rapidity|
875Selector SelectorAbsEtaMin(double absetamin) {
876 return Selector(new SW_QuantityMin<QuantityAbsEta>(absetamin));
877}
878
879// returns a selector for a maximum |pseudo-rapidity|
880Selector SelectorAbsEtaMax(double absetamax) {
881 return Selector(new SW_QuantityMax<QuantityAbsEta>(absetamax));
882}
883
884// returns a selector for a |pseudo-rapidity| range
885Selector SelectorAbsEtaRange(double absetamin, double absetamax) {
886 return Selector(new SW_QuantityRange<QuantityAbsEta>(absetamin, absetamax));
887}
888
889
890//----------------------------------------------------------------------
891/// helper for selecting on azimuthal angle
892///
893/// Note that the bounds have to be specified as min<max
894/// phimin has to be > -2pi
895/// phimax has to be < 4pi
896class SW_PhiRange : public SelectorWorker {
897public:
898 /// detfault ctor (initialises the pt cut)
899 SW_PhiRange(double phimin, double phimax) : _phimin(phimin), _phimax(phimax){
900 assert(_phimin<_phimax);
901 assert(_phimin>-twopi);
902 assert(_phimax<2*twopi);
903
904 _phispan = _phimax - _phimin;
905 }
906
907 /// returns true is the given object passes the selection pt cut
908 virtual bool pass(const PseudoJet & jet) const {
909 double dphi=jet.phi()-_phimin;
910 if (dphi >= twopi) dphi -= twopi;
911 if (dphi < 0) dphi += twopi;
912 return (dphi <= _phispan);
913 }
914
915 /// returns a description of the worker
916 virtual string description() const {
917 ostringstream ostr;
918 ostr << _phimin << " <= phi <= " << _phimax;
919 return ostr.str();
920 }
921
922 virtual bool is_geometric() const { return true;}
923
924protected:
925 double _phimin; // the lower cut
926 double _phimax; // the upper cut
927 double _phispan; // the span of the range
928};
929
930
931// returns a selector for a phi range
932Selector SelectorPhiRange(double phimin, double phimax) {
933 return Selector(new SW_PhiRange(phimin, phimax));
934}
935
936//----------------------------------------------------------------------
937/// helper for selecting on both rapidity and azimuthal angle
938class SW_RapPhiRange : public SW_And{
939public:
940 SW_RapPhiRange(double rapmin, double rapmax, double phimin, double phimax)
941 : SW_And(SelectorRapRange(rapmin, rapmax), SelectorPhiRange(phimin, phimax)){
942 _known_area = ((phimax-phimin > twopi) ? twopi : phimax-phimin) * (rapmax-rapmin);
943 }
944
945 /// if it has a computable area, return it
946 virtual double known_area() const{
947 return _known_area;
948 }
949
950protected:
951 double _known_area;
952};
953
954Selector SelectorRapPhiRange(double rapmin, double rapmax, double phimin, double phimax) {
955 return Selector(new SW_RapPhiRange(rapmin, rapmax, phimin, phimax));
956}
957
958
959//----------------------------------------------------------------------
960/// helper for selecting the n hardest jets
961class SW_NHardest : public SelectorWorker {
962public:
963 /// ctor with specification of the number of objects to keep
964 SW_NHardest(unsigned int n) : _n(n) {};
965
966 /// pass makes no sense here normally the parent selector will throw
967 /// an error but for internal use in the SW, we'll throw one from
968 /// here by security
969 virtual bool pass(const PseudoJet &) const {
970 if (!applies_jet_by_jet())
971 throw Error("Cannot apply this selector worker to an individual jet");
972 return false;
973 }
974
975 /// For each jet that does not pass the cuts, this routine sets the
976 /// pointer to 0.
977 virtual void terminator(vector<const PseudoJet *> & jets) const {
978 // nothing to do if the size is too small
979 if (jets.size() < _n) return;
980
981 // do we want to first chech if things are already ordered before
982 // going through the ordering process? For now, no. Maybe carry
983 // out timing tests at some point to establish the optimal
984 // strategy.
985
986 vector<double> minus_pt2(jets.size());
987 vector<unsigned int> indices(jets.size());
988
989 for (unsigned int i=0; i<jets.size(); i++){
990 indices[i] = i;
991
992 // we need to make sure that the object has not already been
993 // nullified. Note that if we have less than _n jets, this
994 // whole n-hardest selection will not have any effect.
995 minus_pt2[i] = jets[i] ? -jets[i]->perp2() : 0.0;
996 }
997
998 IndexedSortHelper sort_helper(& minus_pt2);
999
1000 partial_sort(indices.begin(), indices.begin()+_n, indices.end(), sort_helper);
1001
1002 for (unsigned int i=_n; i<jets.size(); i++)
1003 jets[indices[i]] = NULL;
1004 }
1005
1006 /// returns true if this can be applied jet by jet
1007 virtual bool applies_jet_by_jet() const {return false;}
1008
1009 /// returns a description of the worker
1010 virtual string description() const {
1011 ostringstream ostr;
1012 ostr << _n << " hardest";
1013 return ostr.str();
1014 }
1015
1016protected:
1017 unsigned int _n;
1018};
1019
1020
1021// returns a selector for the n hardest jets
1022Selector SelectorNHardest(unsigned int n) {
1023 return Selector(new SW_NHardest(n));
1024}
1025
1026
1027
1028//----------------------------------------------------------------------
1029// selector and workers for geometric ranges
1030//----------------------------------------------------------------------
1031
1032//----------------------------------------------------------------------
1033/// a generic class for objects that contain a position
1034class SW_WithReference : public SelectorWorker{
1035public:
1036 /// ctor
1037 SW_WithReference() : _is_initialised(false){};
1038
1039 /// returns true if the worker takes a reference jet
1040 virtual bool takes_reference() const { return true;}
1041
1042 /// sets the reference jet
1043 virtual void set_reference(const PseudoJet &centre){
1044 _is_initialised = true;
1045 _reference = centre;
1046 }
1047
1048protected:
1049 PseudoJet _reference;
1050 bool _is_initialised;
1051};
1052
1053//----------------------------------------------------------------------
1054/// helper for selecting on objects within a distance 'radius' of a reference
1055class SW_Circle : public SW_WithReference {
1056public:
1057 SW_Circle(const double &radius) : _radius2(radius*radius) {}
1058
1059 /// return a copy of the current object
1060 virtual SelectorWorker* copy(){ return new SW_Circle(*this);}
1061
1062 /// returns true if a given object passes the selection criterium
1063 /// this has to be overloaded by derived workers
1064 virtual bool pass(const PseudoJet & jet) const {
1065 // make sure the centre is initialised
1066 if (! _is_initialised)
1067 throw Error("To use a SelectorCircle (or any selector that requires a reference), you first have to call set_reference(...)");
1068
1069 return jet.squared_distance(_reference) <= _radius2;
1070 }
1071
1072 /// returns a description of the worker
1073 virtual string description() const {
1074 ostringstream ostr;
1075 ostr << "distance from the centre <= " << sqrt(_radius2);
1076 return ostr.str();
1077 }
1078
1079 /// returns the rapidity range for which it may return "true"
1080 virtual void get_rapidity_extent(double & rapmin, double & rapmax) const{
1081 // make sure the centre is initialised
1082 if (! _is_initialised)
1083 throw Error("To use a SelectorCircle (or any selector that requires a reference), you first have to call set_reference(...)");
1084
1085 rapmax = _reference.rap()+sqrt(_radius2);
1086 rapmin = _reference.rap()-sqrt(_radius2);
1087 }
1088
1089 virtual bool is_geometric() const { return true;} ///< implies a finite area
1090 virtual bool has_finite_area() const { return true;} ///< regardless of the reference
1091 virtual bool has_known_area() const { return true;} ///< the area is analytically known
1092 virtual double known_area() const {
1093 return pi * _radius2;
1094 }
1095
1096protected:
1097 double _radius2;
1098};
1099
1100
1101// select on objets within a distance 'radius' of a variable location
1102Selector SelectorCircle(const double & radius) {
1103 return Selector(new SW_Circle(radius));
1104}
1105
1106
1107//----------------------------------------------------------------------
1108/// helper for selecting on objects with a distance to a reference
1109/// betwene 'radius_in' and 'radius_out'
1110class SW_Doughnut : public SW_WithReference {
1111public:
1112 SW_Doughnut(const double &radius_in, const double &radius_out)
1113 : _radius_in2(radius_in*radius_in), _radius_out2(radius_out*radius_out) {}
1114
1115 /// return a copy of the current object
1116 virtual SelectorWorker* copy(){ return new SW_Doughnut(*this);}
1117
1118 /// returns true if a given object passes the selection criterium
1119 /// this has to be overloaded by derived workers
1120 virtual bool pass(const PseudoJet & jet) const {
1121 // make sure the centre is initialised
1122 if (! _is_initialised)
1123 throw Error("To use a SelectorDoughnut (or any selector that requires a reference), you first have to call set_reference(...)");
1124
1125 double distance2 = jet.squared_distance(_reference);
1126
1127 return (distance2 <= _radius_out2) && (distance2 >= _radius_in2);
1128 }
1129
1130 /// returns a description of the worker
1131 virtual string description() const {
1132 ostringstream ostr;
1133 ostr << sqrt(_radius_in2) << " <= distance from the centre <= " << sqrt(_radius_out2);
1134 return ostr.str();
1135 }
1136
1137 /// returns the rapidity range for which it may return "true"
1138 virtual void get_rapidity_extent(double & rapmin, double & rapmax) const{
1139 // make sure the centre is initialised
1140 if (! _is_initialised)
1141 throw Error("To use a SelectorDoughnut (or any selector that requires a reference), you first have to call set_reference(...)");
1142
1143 rapmax = _reference.rap()+sqrt(_radius_out2);
1144 rapmin = _reference.rap()-sqrt(_radius_out2);
1145 }
1146
1147 virtual bool is_geometric() const { return true;} ///< implies a finite area
1148 virtual bool has_finite_area() const { return true;} ///< regardless of the reference
1149 virtual bool has_known_area() const { return true;} ///< the area is analytically known
1150 virtual double known_area() const {
1151 return pi * (_radius_out2-_radius_in2);
1152 }
1153
1154protected:
1155 double _radius_in2, _radius_out2;
1156};
1157
1158
1159
1160// select on objets with distance from the centre is between 'radius_in' and 'radius_out'
1161Selector SelectorDoughnut(const double & radius_in, const double & radius_out) {
1162 return Selector(new SW_Doughnut(radius_in, radius_out));
1163}
1164
1165
1166//----------------------------------------------------------------------
1167/// helper for selecting on objects with rapidity within a distance 'delta' of a reference
1168class SW_Strip : public SW_WithReference {
1169public:
1170 SW_Strip(const double &delta) : _delta(delta) {}
1171
1172 /// return a copy of the current object
1173 virtual SelectorWorker* copy(){ return new SW_Strip(*this);}
1174
1175 /// returns true if a given object passes the selection criterium
1176 /// this has to be overloaded by derived workers
1177 virtual bool pass(const PseudoJet & jet) const {
1178 // make sure the centre is initialised
1179 if (! _is_initialised)
1180 throw Error("To use a SelectorStrip (or any selector that requires a reference), you first have to call set_reference(...)");
1181
1182 return abs(jet.rap()-_reference.rap()) <= _delta;
1183 }
1184
1185 /// returns a description of the worker
1186 virtual string description() const {
1187 ostringstream ostr;
1188 ostr << "|rap - rap_reference| <= " << _delta;
1189 return ostr.str();
1190 }
1191
1192 /// returns the rapidity range for which it may return "true"
1193 virtual void get_rapidity_extent(double & rapmin, double & rapmax) const{
1194 // make sure the centre is initialised
1195 if (! _is_initialised)
1196 throw Error("To use a SelectorStrip (or any selector that requires a reference), you first have to call set_reference(...)");
1197
1198 rapmax = _reference.rap()+_delta;
1199 rapmin = _reference.rap()-_delta;
1200 }
1201
1202 virtual bool is_geometric() const { return true;} ///< implies a finite area
1203 virtual bool has_finite_area() const { return true;} ///< regardless of the reference
1204 virtual bool has_known_area() const { return true;} ///< the area is analytically known
1205 virtual double known_area() const {
1206 return twopi * 2 * _delta;
1207 }
1208
1209protected:
1210 double _delta;
1211};
1212
1213
1214// select on objets within a distance 'radius' of a variable location
1215Selector SelectorStrip(const double & half_width) {
1216 return Selector(new SW_Strip(half_width));
1217}
1218
1219
1220//----------------------------------------------------------------------
1221/// helper for selecting on objects with rapidity within a distance
1222/// 'delta_rap' of a reference and phi within a distanve delta_phi of
1223/// a reference
1224class SW_Rectangle : public SW_WithReference {
1225public:
1226 SW_Rectangle(const double &delta_rap, const double &delta_phi)
1227 : _delta_rap(delta_rap), _delta_phi(delta_phi) {}
1228
1229 /// return a copy of the current object
1230 virtual SelectorWorker* copy(){ return new SW_Rectangle(*this);}
1231
1232 /// returns true if a given object passes the selection criterium
1233 /// this has to be overloaded by derived workers
1234 virtual bool pass(const PseudoJet & jet) const {
1235 // make sure the centre is initialised
1236 if (! _is_initialised)
1237 throw Error("To use a SelectorRectangle (or any selector that requires a reference), you first have to call set_reference(...)");
1238
1239 return (abs(jet.rap()-_reference.rap()) <= _delta_rap) && (abs(jet.delta_phi_to(_reference)) <= _delta_phi);
1240 }
1241
1242 /// returns a description of the worker
1243 virtual string description() const {
1244 ostringstream ostr;
1245 ostr << "|rap - rap_reference| <= " << _delta_rap << " && |phi - phi_reference| <= " << _delta_phi ;
1246 return ostr.str();
1247 }
1248
1249 /// returns the rapidity range for which it may return "true"
1250 virtual void get_rapidity_extent(double & rapmin, double & rapmax) const{
1251 // make sure the centre is initialised
1252 if (! _is_initialised)
1253 throw Error("To use a SelectorRectangle (or any selector that requires a reference), you first have to call set_reference(...)");
1254
1255 rapmax = _reference.rap()+_delta_rap;
1256 rapmin = _reference.rap()-_delta_rap;
1257 }
1258
1259 virtual bool is_geometric() const { return true;} ///< implies a finite area
1260 virtual bool has_finite_area() const { return true;} ///< regardless of the reference
1261 virtual bool has_known_area() const { return true;} ///< the area is analytically known
1262 virtual double known_area() const {
1263 return 4 * _delta_rap * _delta_phi;
1264 }
1265
1266protected:
1267 double _delta_rap, _delta_phi;
1268};
1269
1270
1271// select on objets within a distance 'radius' of a variable location
1272Selector SelectorRectangle(const double & half_rap_width, const double & half_phi_width) {
1273 return Selector(new SW_Rectangle(half_rap_width, half_phi_width));
1274}
1275
1276
1277//----------------------------------------------------------------------
1278/// helper for selecting the jets that carry at least a given fraction
1279/// of the reference jet
1280class SW_PtFractionMin : public SW_WithReference {
1281public:
1282 /// ctor with specification of the number of objects to keep
1283 SW_PtFractionMin(double fraction) : _fraction2(fraction*fraction){}
1284
1285 /// return a copy of the current object
1286 virtual SelectorWorker* copy(){ return new SW_PtFractionMin(*this);}
1287
1288 /// return true if the jet carries a large enough fraction of the reference.
1289 /// Throw an error if the reference is not initialised.
1290 virtual bool pass(const PseudoJet & jet) const {
1291 // make sure the centre is initialised
1292 if (! _is_initialised)
1293 throw Error("To use a SelectorPtFractionMin (or any selector that requires a reference), you first have to call set_reference(...)");
1294
1295 // otherwise, just call that method on the jet
1296 return (jet.perp2() >= _fraction2*_reference.perp2());
1297 }
1298
1299 /// returns a description of the worker
1300 virtual string description() const {
1301 ostringstream ostr;
1302 ostr << "pt >= " << sqrt(_fraction2) << "* pt_ref";
1303 return ostr.str();
1304 }
1305
1306protected:
1307 double _fraction2;
1308};
1309
1310
1311// select objects that carry at least a fraction "fraction" of the reference jet
1312// (Note that this selectir takes a reference)
1313Selector SelectorPtFractionMin(double fraction){
1314 return Selector(new SW_PtFractionMin(fraction));
1315}
1316
1317
1318//----------------------------------------------------------------------
1319// additional (mostly helper) selectors
1320//----------------------------------------------------------------------
1321
1322//----------------------------------------------------------------------
1323/// helper for selecting the 0-momentum jets
1324class SW_IsZero : public SelectorWorker {
1325public:
1326 /// ctor
1327 SW_IsZero(){}
1328
1329 /// return true if the jet has zero momentum
1330 virtual bool pass(const PseudoJet & jet) const {
1331 return jet==0;
1332 }
1333
1334 /// rereturns a description of the worker
1335 virtual string description() const { return "zero";}
1336};
1337
1338
1339// select objects with zero momentum
1340Selector SelectorIsZero(){
1341 return Selector(new SW_IsZero());
1342}
1343
1344
1345//----------------------------------------------------------------------
1346#ifndef __FJCORE__
1347/// helper for selecting the pure ghost
1348class SW_IsPureGhost : public SelectorWorker {
1349public:
1350 /// ctor
1351 SW_IsPureGhost(){}
1352
1353 /// return true if the jet is a pure-ghost jet
1354 virtual bool pass(const PseudoJet & jet) const {
1355 // if the jet has no area support then it's certainly not a ghost
1356 if (!jet.has_area()) return false;
1357
1358 // otherwise, just call that method on the jet
1359 return jet.is_pure_ghost();
1360 }
1361
1362 /// rereturns a description of the worker
1363 virtual string description() const { return "pure ghost";}
1364};
1365
1366
1367// select objects that are (or are only made of) ghosts
1368Selector SelectorIsPureGhost(){
1369 return Selector(new SW_IsPureGhost());
1370}
1371
1372//----------------------------------------------------------------------
1373// Selector and workers for obtaining a Selector from an old
1374// RangeDefinition
1375//
1376// This is mostly intended for backward compatibility and is likely to
1377// be removed in a future major release of FastJet
1378//----------------------------------------------------------------------
1379
1380//----------------------------------------------------------------------
1381/// helper for selecting on both rapidity and azimuthal angle
1382class SW_RangeDefinition : public SelectorWorker{
1383public:
1384 /// ctor from a RangeDefinition
1385 SW_RangeDefinition(const RangeDefinition &range) : _range(&range){}
1386
1387 /// transfer the selection creterium to the underlying RangeDefinition
1388 virtual bool pass(const PseudoJet & jet) const {
1389 return _range->is_in_range(jet);
1390 }
1391
1392 /// returns a description of the worker
1393 virtual string description() const {
1394 return _range->description();
1395 }
1396
1397 /// returns the rapidity range for which it may return "true"
1398 virtual void get_rapidity_extent(double & rapmin, double & rapmax) const{
1399 _range->get_rap_limits(rapmin, rapmax);
1400 }
1401
1402 /// check if it has a finite area
1403 virtual bool is_geometric() const { return true;}
1404
1405 /// check if it has an analytically computable area
1406 virtual bool has_known_area() const { return true;}
1407
1408 /// if it has a computable area, return it
1409 virtual double known_area() const{
1410 return _range->area();
1411 }
1412
1413protected:
1414 const RangeDefinition *_range;
1415};
1416
1417
1418// ctor from a RangeDefinition
1419//----------------------------------------------------------------------
1420//
1421// This is provided for backward compatibility and will be removed in
1422// a future major release of FastJet
1423Selector::Selector(const RangeDefinition &range) {
1424 _worker.reset(new SW_RangeDefinition(range));
1425}
1426#endif // __FJCORE__
1427
1428
1429// operators applying directly on a Selector
1430//----------------------------------------------------------------------
1431
1432// operator &=
1433// For 2 Selectors a and b, a &= b is eauivalent to a = a & b;
1434Selector & Selector::operator &=(const Selector & b){
1435 _worker.reset(new SW_And(*this, b));
1436 return *this;
1437}
1438
1439// operator &=
1440// For 2 Selectors a and b, a &= b is eauivalent to a = a & b;
1441Selector & Selector::operator |=(const Selector & b){
1442 _worker.reset(new SW_Or(*this, b));
1443 return *this;
1444}
1445
1446FASTJET_END_NAMESPACE // defined in fastjet/internal/base.hh
Note: See TracBrowser for help on using the repository browser.