34 #include "fastjet/Selector.hh" 
   36 #include "fastjet/GhostedAreaSpec.hh"   
   41 FASTJET_BEGIN_NAMESPACE      
 
   48 std::vector<PseudoJet> Selector::operator()(
const std::vector<PseudoJet> & jets)
 const {
 
   49   std::vector<PseudoJet> result;
 
   54     for (std::vector<PseudoJet>::const_iterator jet = jets.begin(); 
 
   55          jet != jets.end(); jet++) {
 
   56       if (worker_local->
pass(*jet)) result.push_back(*jet);
 
   61     std::vector<const PseudoJet *> jetptrs(jets.size());
 
   62     for (
unsigned i = 0; i < jets.size(); i++) {
 
   63       jetptrs[i] = & jets[i];
 
   66     for (
unsigned i = 0; i < jetptrs.size(); i++) {
 
   67       if (jetptrs[i]) result.push_back(jets[i]);
 
   76 unsigned int Selector::count(
const std::vector<PseudoJet> & jets)
 const {
 
   82     for (
unsigned i = 0; i < jets.size(); i++) {
 
   83       if (worker_local->
pass(jets[i])) n++;
 
   86     std::vector<const PseudoJet *> jetptrs(jets.size());
 
   87     for (
unsigned i = 0; i < jets.size(); i++) {
 
   88       jetptrs[i] = & jets[i];
 
   91     for (
unsigned i = 0; i < jetptrs.size(); i++) {
 
  101 PseudoJet Selector::sum(
const std::vector<PseudoJet> & jets)
 const {
 
  107     for (
unsigned i = 0; i < jets.size(); i++) {
 
  108       if (worker_local->
pass(jets[i])) this_sum += jets[i];
 
  111     std::vector<const PseudoJet *> jetptrs(jets.size());
 
  112     for (
unsigned i = 0; i < jets.size(); i++) {
 
  113       jetptrs[i] = & jets[i];
 
  116     for (
unsigned i = 0; i < jetptrs.size(); i++) {
 
  117       if (jetptrs[i]) this_sum += jets[i];
 
  126 double Selector::scalar_pt_sum(
const std::vector<PseudoJet> & jets)
 const {
 
  127   double this_sum = 0.0;
 
  132     for (
unsigned i = 0; i < jets.size(); i++) {
 
  133       if (worker_local->
pass(jets[i])) this_sum += jets[i].pt();
 
  136     std::vector<const PseudoJet *> jetptrs(jets.size());
 
  137     for (
unsigned i = 0; i < jets.size(); i++) {
 
  138       jetptrs[i] = & jets[i];
 
  141     for (
unsigned i = 0; i < jetptrs.size(); i++) {
 
  142       if (jetptrs[i]) this_sum += jets[i].pt();
 
  153 void Selector::sift(
const std::vector<PseudoJet> & jets,
 
  154                     std::vector<PseudoJet> & jets_that_pass,
 
  155                     std::vector<PseudoJet> & jets_that_fail
 
  159   jets_that_pass.clear();
 
  160   jets_that_fail.clear();
 
  164     for (
unsigned i = 0; i < jets.size(); i++) {
 
  165       if (worker_local->
pass(jets[i])) {
 
  166         jets_that_pass.push_back(jets[i]);
 
  168         jets_that_fail.push_back(jets[i]);
 
  172     std::vector<const PseudoJet *> jetptrs(jets.size());
 
  173     for (
unsigned i = 0; i < jets.size(); i++) {
 
  174       jetptrs[i] = & jets[i];
 
  177     for (
unsigned i = 0; i < jetptrs.size(); i++) {
 
  179         jets_that_pass.push_back(jets[i]);
 
  181         jets_that_fail.push_back(jets[i]);
 
  189 double Selector::area()
 const{
 
  190   return area(gas::def_ghost_area);
 
  194 double Selector::area(
double ghost_area)
 const{
 
  198   if (_worker->has_known_area()) 
return _worker->known_area();
 
  201   double rapmin, rapmax;
 
  202   get_rapidity_extent(rapmin, rapmax);
 
  204   std::vector<PseudoJet> ghosts;
 
  208   return ghost_spec.ghost_area() * ((*this)(ghosts)).size();
 
  217 bool SelectorWorker::has_finite_area()
 const { 
 
  218   if (! is_geometric()) 
return false;
 
  219   double rapmin, rapmax;
 
  220   get_rapidity_extent(rapmin, rapmax);
 
  221   return (rapmax != std::numeric_limits<double>::infinity())
 
  222     &&  (-rapmin != std::numeric_limits<double>::infinity());
 
  239   virtual bool pass(
const PseudoJet &)
 const {
 
  245   virtual void terminator(vector<const PseudoJet *> &)
 const {
 
  251   virtual string description()
 const { 
return "Identity";}
 
  254   virtual bool is_geometric()
 const { 
return true;}
 
  259 Selector SelectorIdentity() {
 
  260   return Selector(
new SW_Identity);
 
  270 class SW_Not : 
public SelectorWorker {
 
  273   SW_Not(
const Selector & s) : _s(s) {}
 
  276   virtual SelectorWorker* copy(){ 
return new SW_Not(*
this);}
 
  280   virtual bool pass(
const PseudoJet & jet)
 const {
 
  282     if (!applies_jet_by_jet())
 
  283       throw Error(
"Cannot apply this selector worker to an individual jet");
 
  285     return ! _s.pass(jet);
 
  289   virtual bool applies_jet_by_jet()
 const {
return _s.applies_jet_by_jet();}
 
  292   virtual void terminator(vector<const PseudoJet *> & jets)
 const {
 
  295     if (applies_jet_by_jet()){
 
  296       SelectorWorker::terminator(jets);
 
  301     vector<const PseudoJet *> s_jets = jets;
 
  302     _s.worker()->terminator(s_jets);
 
  306     for (
unsigned int i=0; i<s_jets.size(); i++){
 
  307       if (s_jets[i]) jets[i] = NULL;
 
  312   virtual string description()
 const {
 
  314     ostr << 
"!(" << _s.description() << 
")";
 
  319   virtual bool is_geometric()
 const { 
return _s.is_geometric();}
 
  322   virtual bool takes_reference()
 const { 
return _s.takes_reference();}
 
  325   virtual void set_reference(
const PseudoJet &ref) { _s.set_reference(ref);}
 
  340 class SW_BinaryOperator: 
public SelectorWorker {
 
  343   SW_BinaryOperator(
const Selector & s1, 
const Selector & s2) : _s1(s1), _s2(s2) {
 
  347     _applies_jet_by_jet = _s1.applies_jet_by_jet() && _s2.applies_jet_by_jet();
 
  350     _takes_reference = _s1.takes_reference() || _s2.takes_reference();
 
  353     _is_geometric = _s1.is_geometric() && _s2.is_geometric();
 
  357   virtual bool applies_jet_by_jet()
 const {
return _applies_jet_by_jet;}
 
  360   virtual bool takes_reference()
 const{ 
 
  361     return _takes_reference;
 
  365   virtual void set_reference(
const PseudoJet ¢re){
 
  366     _s1.set_reference(centre);
 
  367     _s2.set_reference(centre);
 
  371   virtual bool is_geometric()
 const { 
return _is_geometric;} 
 
  375   bool _applies_jet_by_jet;
 
  376   bool _takes_reference;
 
  384 class SW_And: 
public SW_BinaryOperator {
 
  387   SW_And(
const Selector & s1, 
const Selector & s2) : SW_BinaryOperator(s1,s2){}
 
  390   virtual SelectorWorker* copy(){ 
return new SW_And(*
this);}
 
  394   virtual bool pass(
const PseudoJet & jet)
 const {
 
  396     if (!applies_jet_by_jet())
 
  397       throw Error(
"Cannot apply this selector worker to an individual jet");
 
  399     return _s1.pass(jet) && _s2.pass(jet);
 
  403   virtual void terminator(vector<const PseudoJet *> & jets)
 const {
 
  406     if (applies_jet_by_jet()){
 
  407       SelectorWorker::terminator(jets);
 
  412     vector<const PseudoJet *> s1_jets = jets;
 
  413     _s1.worker()->terminator(s1_jets);
 
  416     _s2.worker()->terminator(jets);
 
  419     for (
unsigned int i=0; i<jets.size(); i++){
 
  420       if (! s1_jets[i]) jets[i] = NULL;
 
  425   virtual void get_rapidity_extent(
double & rapmin, 
double & rapmax)
 const {
 
  426     double s1min, s1max, s2min, s2max;
 
  427     _s1.get_rapidity_extent(s1min, s1max);
 
  428     _s2.get_rapidity_extent(s2min, s2max);
 
  429     rapmax = min(s1max, s2max);
 
  430     rapmin = max(s1min, s2min);
 
  434   virtual string description()
 const {
 
  436     ostr << 
"(" << _s1.description() << 
" && " << _s2.description() << 
")";
 
  451 class SW_Or: 
public SW_BinaryOperator {
 
  454   SW_Or(
const Selector & s1, 
const Selector & s2) : SW_BinaryOperator(s1,s2) {}
 
  457   virtual SelectorWorker* copy(){ 
return new SW_Or(*
this);}
 
  461   virtual bool pass(
const PseudoJet & jet)
 const {
 
  463     if (!applies_jet_by_jet())
 
  464       throw Error(
"Cannot apply this selector worker to an individual jet");
 
  466     return _s1.pass(jet) || _s2.pass(jet);
 
  470   virtual bool applies_jet_by_jet()
 const {
 
  474     return _s1.applies_jet_by_jet() && _s2.applies_jet_by_jet();
 
  478   virtual void terminator(vector<const PseudoJet *> & jets)
 const {
 
  481     if (applies_jet_by_jet()){
 
  482       SelectorWorker::terminator(jets);
 
  487     vector<const PseudoJet *> s1_jets = jets;
 
  488     _s1.worker()->terminator(s1_jets);
 
  491     _s2.worker()->terminator(jets);
 
  495     for (
unsigned int i=0; i<jets.size(); i++){
 
  496       if (s1_jets[i]) jets[i] = s1_jets[i];
 
  501   virtual string description()
 const {
 
  503     ostr << 
"(" << _s1.description() << 
" || " << _s2.description() << 
")";
 
  508   virtual void get_rapidity_extent(
double & rapmin, 
double & rapmax)
 const {
 
  509     double s1min, s1max, s2min, s2max;
 
  510     _s1.get_rapidity_extent(s1min, s1max);
 
  511     _s2.get_rapidity_extent(s2min, s2max);
 
  512     rapmax = max(s1max, s2max);
 
  513     rapmin = min(s1min, s2min);
 
  525 class SW_Mult: 
public SW_And {
 
  528   SW_Mult(
const Selector & s1, 
const Selector & s2) : SW_And(s1,s2) {}
 
  531   virtual SelectorWorker* copy(){ 
return new SW_Mult(*
this);}
 
  534   virtual void terminator(vector<const PseudoJet *> & jets)
 const {
 
  537     if (applies_jet_by_jet()){
 
  538       SelectorWorker::terminator(jets);
 
  543     _s2.worker()->terminator(jets);
 
  546     _s1.worker()->terminator(jets);
 
  550   virtual string description()
 const {
 
  552     ostr << 
"(" << _s1.description() << 
" * " << _s2.description() << 
")";
 
  560   return Selector(
new SW_Mult(s1,s2));
 
  575   QuantityBase(
double q) : _q(q){}
 
  576   virtual ~QuantityBase(){}
 
  577   virtual double operator()(
const PseudoJet & jet ) 
const =0;
 
  578   virtual string description() 
const =0;
 
  579   virtual bool is_geometric()
 const { 
return false;}
 
  580   virtual double comparison_value()
 const {
return _q;}
 
  581   virtual double description_value()
 const {
return comparison_value();}
 
  587 class QuantitySquareBase : 
public QuantityBase{
 
  589   QuantitySquareBase(
double sqrtq) : QuantityBase(sqrtq*sqrtq), _sqrtq(sqrtq){}
 
  590   virtual double description_value()
 const {
return _sqrtq;}
 
  596 template<
typename QuantityType>
 
  597 class SW_QuantityMin : 
public SelectorWorker{
 
  600   SW_QuantityMin(
double qmin) : _qmin(qmin) {}
 
  603   virtual bool pass(
const PseudoJet & jet)
 const {
return _qmin(jet) >= _qmin.comparison_value();}
 
  606   virtual string description()
 const {
 
  608     ostr << _qmin.description() << 
" >= " << _qmin.description_value();
 
  612   virtual bool is_geometric()
 const { 
return _qmin.is_geometric();}
 
  620 template<
typename QuantityType>
 
  621 class SW_QuantityMax : 
public SelectorWorker {
 
  624   SW_QuantityMax(
double qmax) : _qmax(qmax) {}
 
  627   virtual bool pass(
const PseudoJet & jet)
 const {
return _qmax(jet) <= _qmax.comparison_value();}
 
  630   virtual string description()
 const {
 
  632     ostr << _qmax.description() << 
" <= " << _qmax.description_value();
 
  636   virtual bool is_geometric()
 const { 
return _qmax.is_geometric();}
 
  644 template<
typename QuantityType>
 
  645 class SW_QuantityRange : 
public SelectorWorker {
 
  648   SW_QuantityRange(
double qmin, 
double qmax) : _qmin(qmin), _qmax(qmax) {}
 
  651   virtual bool pass(
const PseudoJet & jet)
 const {
 
  652     double q = _qmin(jet); 
 
  653     return (q >= _qmin.comparison_value()) && (q <= _qmax.comparison_value());
 
  657   virtual string description()
 const {
 
  659     ostr << _qmin.description_value() << 
" <= " << _qmin.description() << 
" <= " << _qmax.description_value();
 
  663   virtual bool is_geometric()
 const { 
return _qmin.is_geometric();}
 
  673 class QuantityPt2 : 
public QuantitySquareBase{
 
  675   QuantityPt2(
double pt) : QuantitySquareBase(pt){}
 
  676   virtual double operator()(
const PseudoJet & jet )
 const { 
return jet.perp2();}
 
  677   virtual string description()
 const {
return "pt";}
 
  682   return Selector(
new SW_QuantityMin<QuantityPt2>(ptmin));
 
  687   return Selector(
new SW_QuantityMax<QuantityPt2>(ptmax));
 
  692   return Selector(
new SW_QuantityRange<QuantityPt2>(ptmin, ptmax));
 
  698 class QuantityEt2 : 
public QuantitySquareBase{
 
  700   QuantityEt2(
double Et) : QuantitySquareBase(Et){}
 
  701   virtual double operator()(
const PseudoJet & jet )
 const { 
return jet.Et2();}
 
  702   virtual string description()
 const {
return "Et";}
 
  707   return Selector(
new SW_QuantityMin<QuantityEt2>(Etmin));
 
  712   return Selector(
new SW_QuantityMax<QuantityEt2>(Etmax));
 
  717   return Selector(
new SW_QuantityRange<QuantityEt2>(Etmin, Etmax));
 
  723 class QuantityE : 
public QuantityBase{
 
  725   QuantityE(
double E) : QuantityBase(E){}
 
  726   virtual double operator()(
const PseudoJet & jet )
 const { 
return jet.E();}
 
  727   virtual string description()
 const {
return "E";}
 
  732   return Selector(
new SW_QuantityMin<QuantityE>(Emin));
 
  737   return Selector(
new SW_QuantityMax<QuantityE>(Emax));
 
  742   return Selector(
new SW_QuantityRange<QuantityE>(Emin, Emax));
 
  748 class QuantityM2 : 
public QuantitySquareBase{
 
  750   QuantityM2(
double m) : QuantitySquareBase(m){}
 
  751   virtual double operator()(
const PseudoJet & jet )
 const { 
return jet.m2();}
 
  752   virtual string description()
 const {
return "mass";}
 
  757   return Selector(
new SW_QuantityMin<QuantityM2>(mmin));
 
  762   return Selector(
new SW_QuantityMax<QuantityM2>(mmax));
 
  767   return Selector(
new SW_QuantityRange<QuantityM2>(mmin, mmax));
 
  774 class QuantityRap : 
public QuantityBase{
 
  776   QuantityRap(
double rap) : QuantityBase(rap){}
 
  777   virtual double operator()(
const PseudoJet & jet )
 const { 
return jet.rap();}
 
  778   virtual string description()
 const {
return "rap";}
 
  779   virtual bool is_geometric()
 const { 
return true;}
 
  784 class SW_RapMin : 
public SW_QuantityMin<QuantityRap>{
 
  786   SW_RapMin(
double rapmin) : SW_QuantityMin<QuantityRap>(rapmin){}
 
  787   virtual void get_rapidity_extent(
double &rapmin, 
double & rapmax)
 const{
 
  788     rapmax = std::numeric_limits<double>::max();     
 
  789     rapmin = _qmin.comparison_value();
 
  794 class SW_RapMax : 
public SW_QuantityMax<QuantityRap>{
 
  796   SW_RapMax(
double rapmax) : SW_QuantityMax<QuantityRap>(rapmax){}
 
  797   virtual void get_rapidity_extent(
double &rapmin, 
double & rapmax)
 const{
 
  798     rapmax = _qmax.comparison_value(); 
 
  799     rapmin = -std::numeric_limits<double>::max();
 
  804 class SW_RapRange : 
public SW_QuantityRange<QuantityRap>{
 
  806   SW_RapRange(
double rapmin, 
double rapmax) : SW_QuantityRange<QuantityRap>(rapmin, rapmax){
 
  807     assert(rapmin<=rapmax);
 
  809   virtual void get_rapidity_extent(
double &rapmin, 
double & rapmax)
 const{
 
  810     rapmax = _qmax.comparison_value();      
 
  811     rapmin = _qmin.comparison_value(); 
 
  813   virtual bool has_known_area()
 const { 
return true;} 
 
  814   virtual double known_area()
 const { 
 
  815     return twopi * (_qmax.comparison_value()-_qmin.comparison_value());
 
  821   return Selector(
new SW_RapMin(rapmin));
 
  826   return Selector(
new SW_RapMax(rapmax));
 
  831   return Selector(
new SW_RapRange(rapmin, rapmax));
 
  837 class QuantityAbsRap : 
public QuantityBase{
 
  839   QuantityAbsRap(
double absrap) : QuantityBase(absrap){}
 
  840   virtual double operator()(
const PseudoJet & jet )
 const { 
return abs(jet.rap());}
 
  841   virtual string description()
 const {
return "|rap|";}
 
  842   virtual bool is_geometric()
 const { 
return true;}
 
  847 class SW_AbsRapMax : 
public SW_QuantityMax<QuantityAbsRap>{
 
  849   SW_AbsRapMax(
double absrapmax) : SW_QuantityMax<QuantityAbsRap>(absrapmax){}
 
  850   virtual void get_rapidity_extent(
double &rapmin, 
double & rapmax)
 const{
 
  851     rapmax =  _qmax.comparison_value(); 
 
  852     rapmin = -_qmax.comparison_value();
 
  854   virtual bool has_known_area()
 const { 
return true;}   
 
  855   virtual double known_area()
 const { 
 
  856     return twopi * 2 * _qmax.comparison_value();
 
  861 class SW_AbsRapRange : 
public SW_QuantityRange<QuantityAbsRap>{
 
  863   SW_AbsRapRange(
double absrapmin, 
double absrapmax) : SW_QuantityRange<QuantityAbsRap>(absrapmin, absrapmax){}
 
  864   virtual void get_rapidity_extent(
double &rapmin, 
double & rapmax)
 const{
 
  865     rapmax =  _qmax.comparison_value(); 
 
  866     rapmin = -_qmax.comparison_value();
 
  868   virtual bool has_known_area()
 const { 
return true;} 
 
  869   virtual double known_area()
 const { 
 
  870     return twopi * 2 * (_qmax.comparison_value()-max(_qmin.comparison_value(),0.0)); 
 
  876   return Selector(
new SW_QuantityMin<QuantityAbsRap>(absrapmin));
 
  881   return Selector(
new SW_AbsRapMax(absrapmax));
 
  886   return Selector(
new SW_AbsRapRange(rapmin, rapmax));
 
  892 class QuantityEta : 
public QuantityBase{
 
  894   QuantityEta(
double eta) : QuantityBase(eta){}
 
  895   virtual double operator()(
const PseudoJet & jet )
 const { 
return jet.eta();}
 
  896   virtual string description()
 const {
return "eta";}
 
  902   return Selector(
new SW_QuantityMin<QuantityEta>(etamin));
 
  907   return Selector(
new SW_QuantityMax<QuantityEta>(etamax));
 
  912   return Selector(
new SW_QuantityRange<QuantityEta>(etamin, etamax));
 
  918 class QuantityAbsEta : 
public QuantityBase{
 
  920   QuantityAbsEta(
double abseta) : QuantityBase(abseta){}
 
  921   virtual double operator()(
const PseudoJet & jet )
 const { 
return abs(jet.eta());}
 
  922   virtual string description()
 const {
return "|eta|";}
 
  923   virtual bool is_geometric()
 const { 
return true;}
 
  928   return Selector(
new SW_QuantityMin<QuantityAbsEta>(absetamin));
 
  933   return Selector(
new SW_QuantityMax<QuantityAbsEta>(absetamax));
 
  938   return Selector(
new SW_QuantityRange<QuantityAbsEta>(absetamin, absetamax));
 
  948 class SW_PhiRange : 
public SelectorWorker {
 
  951   SW_PhiRange(
double phimin, 
double phimax) : _phimin(phimin), _phimax(phimax){
 
  952     assert(_phimin<_phimax);
 
  953     assert(_phimin>-twopi);
 
  954     assert(_phimax<2*twopi);
 
  956     _phispan = _phimax - _phimin;
 
  960   virtual bool pass(
const PseudoJet & jet)
 const {
 
  961     double dphi=jet.phi()-_phimin;
 
  962     if (dphi >= twopi) dphi -= twopi;
 
  963     if (dphi < 0)      dphi += twopi;
 
  964     return (dphi <= _phispan);
 
  968   virtual string description()
 const {
 
  970     ostr << _phimin << 
" <= phi <= " << _phimax;
 
  974   virtual bool is_geometric()
 const { 
return true;}
 
  985   return Selector(
new SW_PhiRange(phimin, phimax));
 
  990 class SW_RapPhiRange : 
public SW_And{
 
  992   SW_RapPhiRange(
double rapmin, 
double rapmax, 
double phimin, 
double phimax)
 
  994     _known_area = ((phimax-phimin > twopi) ? twopi : phimax-phimin) * (rapmax-rapmin);
 
  998   virtual double known_area()
 const{
 
 1007   return Selector(
new SW_RapPhiRange(rapmin, rapmax, phimin, phimax));
 
 1013 class SW_NHardest : 
public SelectorWorker {
 
 1016   SW_NHardest(
unsigned int n) : _n(n) {};
 
 1021   virtual bool pass(
const PseudoJet &)
 const {
 
 1022     if (!applies_jet_by_jet())
 
 1023       throw Error(
"Cannot apply this selector worker to an individual jet");
 
 1029   virtual void terminator(vector<const PseudoJet *> & jets)
 const {
 
 1031     if (jets.size() < _n) 
return;
 
 1038     vector<double> minus_pt2(jets.size());
 
 1039     vector<unsigned int> indices(jets.size());
 
 1041     for (
unsigned int i=0; i<jets.size(); i++){
 
 1047       minus_pt2[i] = jets[i] ? -jets[i]->perp2() : 0.0;
 
 1050     IndexedSortHelper sort_helper(& minus_pt2);
 
 1052     partial_sort(indices.begin(), indices.begin()+_n, indices.end(), sort_helper);
 
 1054     for (
unsigned int i=_n; i<jets.size(); i++)
 
 1055       jets[indices[i]] = NULL;
 
 1059   virtual bool applies_jet_by_jet()
 const {
return false;}
 
 1062   virtual string description()
 const {
 
 1064     ostr << _n << 
" hardest";
 
 1075   return Selector(
new SW_NHardest(n));
 
 1086 class SW_WithReference : 
public SelectorWorker{
 
 1089   SW_WithReference() : _is_initialised(false){};
 
 1092   virtual bool takes_reference()
 const { 
return true;}
 
 1095   virtual void set_reference(
const PseudoJet ¢re){
 
 1096     _is_initialised = 
true;
 
 1097     _reference = centre;
 
 1101   PseudoJet _reference;
 
 1102   bool _is_initialised;
 
 1107 class SW_Circle : 
public SW_WithReference {
 
 1109   SW_Circle(
const double radius) : _radius2(radius*radius) {}
 
 1112   virtual SelectorWorker* copy(){ 
return new SW_Circle(*
this);}
 
 1116   virtual bool pass(
const PseudoJet & jet)
 const {
 
 1118     if (! _is_initialised)
 
 1119       throw Error(
"To use a SelectorCircle (or any selector that requires a reference), you first have to call set_reference(...)");
 
 1121     return jet.squared_distance(_reference) <= _radius2;
 
 1125   virtual string description()
 const {
 
 1127     ostr << 
"distance from the centre <= " << sqrt(_radius2);
 
 1132   virtual void get_rapidity_extent(
double & rapmin, 
double & rapmax)
 const{
 
 1134     if (! _is_initialised)
 
 1135       throw Error(
"To use a SelectorCircle (or any selector that requires a reference), you first have to call set_reference(...)");
 
 1137     rapmax = _reference.rap()+sqrt(_radius2);
 
 1138     rapmin = _reference.rap()-sqrt(_radius2);
 
 1141   virtual bool is_geometric()
 const { 
return true;}    
 
 1142   virtual bool has_finite_area()
 const { 
return true;} 
 
 1143   virtual bool has_known_area()
 const { 
return true;}  
 
 1144   virtual double known_area()
 const { 
 
 1145     return pi * _radius2;
 
 1155   return Selector(
new SW_Circle(radius));
 
 1162 class SW_Doughnut : 
public SW_WithReference {
 
 1164   SW_Doughnut(
const double radius_in, 
const double radius_out)
 
 1165     : _radius_in2(radius_in*radius_in), _radius_out2(radius_out*radius_out) {}
 
 1168   virtual SelectorWorker* copy(){ 
return new SW_Doughnut(*
this);}
 
 1172   virtual bool pass(
const PseudoJet & jet)
 const {
 
 1174     if (! _is_initialised)
 
 1175       throw Error(
"To use a SelectorDoughnut (or any selector that requires a reference), you first have to call set_reference(...)");
 
 1177     double distance2 = jet.squared_distance(_reference);
 
 1179     return (distance2 <= _radius_out2) && (distance2 >= _radius_in2);
 
 1183   virtual string description()
 const {
 
 1185     ostr << sqrt(_radius_in2) << 
" <= distance from the centre <= " << sqrt(_radius_out2);
 
 1190   virtual void get_rapidity_extent(
double & rapmin, 
double & rapmax)
 const{
 
 1192     if (! _is_initialised)
 
 1193       throw Error(
"To use a SelectorDoughnut (or any selector that requires a reference), you first have to call set_reference(...)");
 
 1195     rapmax = _reference.rap()+sqrt(_radius_out2);
 
 1196     rapmin = _reference.rap()-sqrt(_radius_out2);
 
 1199   virtual bool is_geometric()
 const { 
return true;}    
 
 1200   virtual bool has_finite_area()
 const { 
return true;} 
 
 1201   virtual bool has_known_area()
 const { 
return true;}  
 
 1202   virtual double known_area()
 const { 
 
 1203     return pi * (_radius_out2-_radius_in2);
 
 1207   double _radius_in2, _radius_out2;
 
 1214   return Selector(
new SW_Doughnut(radius_in, radius_out));
 
 1220 class SW_Strip : 
public SW_WithReference {
 
 1222   SW_Strip(
const double delta) : _delta(delta) {}
 
 1225   virtual SelectorWorker* copy(){ 
return new SW_Strip(*
this);}
 
 1229   virtual bool pass(
const PseudoJet & jet)
 const {
 
 1231     if (! _is_initialised)
 
 1232       throw Error(
"To use a SelectorStrip (or any selector that requires a reference), you first have to call set_reference(...)");
 
 1234     return abs(jet.rap()-_reference.rap()) <= _delta;
 
 1238   virtual string description()
 const {
 
 1240     ostr << 
"|rap - rap_reference| <= " << _delta;
 
 1245   virtual void get_rapidity_extent(
double & rapmin, 
double & rapmax)
 const{
 
 1247     if (! _is_initialised)
 
 1248       throw Error(
"To use a SelectorStrip (or any selector that requires a reference), you first have to call set_reference(...)");
 
 1250     rapmax = _reference.rap()+_delta;
 
 1251     rapmin = _reference.rap()-_delta;
 
 1254   virtual bool is_geometric()
 const { 
return true;}    
 
 1255   virtual bool has_finite_area()
 const { 
return true;} 
 
 1256   virtual bool has_known_area()
 const { 
return true;}  
 
 1257   virtual double known_area()
 const { 
 
 1258     return twopi * 2 * _delta;
 
 1268   return Selector(
new SW_Strip(half_width));
 
 1276 class SW_Rectangle : 
public SW_WithReference {
 
 1278   SW_Rectangle(
const double delta_rap, 
const double delta_phi)
 
 1279     : _delta_rap(delta_rap),  _delta_phi(delta_phi) {}
 
 1282   virtual SelectorWorker* copy(){ 
return new SW_Rectangle(*
this);}
 
 1286   virtual bool pass(
const PseudoJet & jet)
 const {
 
 1288     if (! _is_initialised)
 
 1289       throw Error(
"To use a SelectorRectangle (or any selector that requires a reference), you first have to call set_reference(...)");
 
 1291     return (abs(jet.rap()-_reference.rap()) <= _delta_rap) && (abs(jet.delta_phi_to(_reference)) <= _delta_phi);
 
 1295   virtual string description()
 const {
 
 1297     ostr << 
"|rap - rap_reference| <= " << _delta_rap << 
" && |phi - phi_reference| <= " << _delta_phi ;
 
 1302   virtual void get_rapidity_extent(
double & rapmin, 
double & rapmax)
 const{
 
 1304     if (! _is_initialised)
 
 1305       throw Error(
"To use a SelectorRectangle (or any selector that requires a reference), you first have to call set_reference(...)");
 
 1307     rapmax = _reference.rap()+_delta_rap;
 
 1308     rapmin = _reference.rap()-_delta_rap;
 
 1311   virtual bool is_geometric()
 const { 
return true;}    
 
 1312   virtual bool has_finite_area()
 const { 
return true;} 
 
 1313   virtual bool has_known_area()
 const { 
return true;}  
 
 1314   virtual double known_area()
 const { 
 
 1315     return 4 * _delta_rap * _delta_phi;
 
 1319   double _delta_rap, _delta_phi;
 
 1325   return Selector(
new SW_Rectangle(half_rap_width, half_phi_width));
 
 1332 class SW_PtFractionMin : 
public SW_WithReference {
 
 1335   SW_PtFractionMin(
double fraction) : _fraction2(fraction*fraction){}
 
 1338   virtual SelectorWorker* copy(){ 
return new SW_PtFractionMin(*
this);}
 
 1342   virtual bool pass(
const PseudoJet & jet)
 const {
 
 1344     if (! _is_initialised)
 
 1345       throw Error(
"To use a SelectorPtFractionMin (or any selector that requires a reference), you first have to call set_reference(...)");
 
 1348     return (jet.perp2() >= _fraction2*_reference.perp2());
 
 1352   virtual string description()
 const {
 
 1354     ostr << 
"pt >= " << sqrt(_fraction2) << 
"* pt_ref";
 
 1366   return Selector(
new SW_PtFractionMin(fraction));
 
 1376 class SW_IsZero : 
public SelectorWorker {
 
 1382   virtual bool pass(
const PseudoJet & jet)
 const {
 
 1387   virtual string description()
 const { 
return "zero";}
 
 1400 class SW_IsPureGhost : 
public SelectorWorker {
 
 1406   virtual bool pass(
const PseudoJet & jet)
 const {
 
 1408     if (!jet.has_area()) 
return false;
 
 1411     return jet.is_pure_ghost();
 
 1415   virtual string description()
 const { 
return "pure ghost";}
 
 1421   return Selector(
new SW_IsPureGhost());
 
 1434 class SW_RangeDefinition : 
public SelectorWorker{
 
 1437   SW_RangeDefinition(
const RangeDefinition &range) : _range(&range){}
 
 1440   virtual bool pass(
const PseudoJet & jet)
 const {
 
 1441     return _range->is_in_range(jet);
 
 1445   virtual string description()
 const {
 
 1446     return _range->description();
 
 1450   virtual void get_rapidity_extent(
double & rapmin, 
double & rapmax)
 const{
 
 1451     _range->get_rap_limits(rapmin, rapmax);
 
 1455   virtual bool is_geometric()
 const { 
return true;}
 
 1458   virtual bool has_known_area()
 const { 
return true;}
 
 1461   virtual double known_area()
 const{
 
 1462     return _range->area();
 
 1466   const RangeDefinition *_range;
 
 1476   _worker.reset(
new SW_RangeDefinition(range));
 
 1478 #endif  // __FJCORE__ 
 1487   _worker.reset(
new SW_And(*
this, b));
 
 1494   _worker.reset(
new SW_Or(*
this, b));
 
 1498 FASTJET_END_NAMESPACE