FastJet 3.0beta1
11-filter.cc
Go to the documentation of this file.
00001 //----------------------------------------------------------------------
00002 /// \file
00003 /// \page Example11 11 - use of filtering
00004 ///
00005 /// fastjet example program to illustrate the use of the fastjet::Filter class
00006 ///
00007 /// We apply different filter examples to either the hardest jet of the given event, 
00008 /// or to the composition of the two hardest jets: 
00009 ///   - two examples of a filter keeping a fixed number of subjets (as in arXiv:0802.2470)
00010 ///   - a "trimmer" i.e. a filter keeping subjets carrying at least a given 
00011 ///     fraction of the pt of the jet (arXiv:0912.1342).
00012 ///   - two examples of filter in combination with background subtraction
00013 ///
00014 /// run it with    : ./11-filter < data/single-event.dat
00015 ///
00016 /// Source code: 11-filter.cc
00017 //----------------------------------------------------------------------
00018 
00019 #include <fastjet/PseudoJet.hh>
00020 #include <fastjet/ClusterSequence.hh>
00021 #include <fastjet/Selector.hh>
00022 #include <iostream>
00023 #include "fastjet/tools/Filter.hh"
00024 // the following includes are only needed when combining filtering with subtraction
00025 #include "fastjet/tools/GridMedianBackgroundEstimator.hh"
00026 #include "fastjet/ClusterSequenceArea.hh"
00027 #include "fastjet/tools/Subtractor.hh"
00028 
00029 #include <cstdio>   // needed for io
00030 
00031 using namespace fastjet;
00032 using namespace std;
00033 
00034 // a function returning
00035 //   min(Rmax, deltaR_factor * deltaR(j1,j2))
00036 // where j1 and j2 are the 2 subjets of j
00037 // if the jet does not have 2 exactly pieces, Rmax is used.
00038 class DynamicRfilt : public FunctionOfPseudoJet<double>{
00039 public:
00040   // default ctor 
00041   DynamicRfilt(double Rmax, double deltaR_factor) : _Rmax(Rmax), _deltaR_factor(deltaR_factor){}
00042 
00043   // action of the function
00044   double result(const PseudoJet &j) const{
00045     if (! j.has_pieces()) return _Rmax;
00046 
00047     vector<PseudoJet> pieces = j.pieces();
00048     if (! pieces.size()==2) return _Rmax;
00049 
00050     double deltaR = pieces[0].delta_R(pieces[1]);
00051     return min(_Rmax, _deltaR_factor * deltaR);
00052   }
00053 
00054 private:
00055   double _Rmax, _deltaR_factor;
00056 };
00057 
00058 /// an example program showing how to use Filter in FastJet
00059 int main (int argc, char ** argv) {
00060   // read in input particles
00061   //----------------------------------------------------------
00062   vector<PseudoJet> input_particles;
00063   
00064   double px, py , pz, E;
00065   while (cin >> px >> py >> pz >> E) {
00066     // create a fastjet::PseudoJet with these components and put it onto
00067     // back of the input_particles vector
00068     input_particles.push_back(PseudoJet(px,py,pz,E)); 
00069   }
00070  
00071   // get the resulting jets ordered in pt
00072   //----------------------------------------------------------
00073   JetDefinition jet_def(cambridge_algorithm, 1.2);
00074   // the use of a ClusterSequenceArea (instead of a plain ClusterSequence)
00075   // is only needed because we will later combine filtering with area-based
00076   // subtraction
00077   ClusterSequenceArea clust_seq(input_particles, jet_def, AreaDefinition(active_area_explicit_ghosts));
00078   vector<fastjet::PseudoJet> inclusive_jets = sorted_by_pt(clust_seq.inclusive_jets(5.0));
00079 
00080   // label the columns
00081   printf("%5s %15s %15s %15s\n","jet #", "rapidity", "phi", "pt");
00082  
00083   // print out the details for each jet
00084   for (unsigned int i = 0; i < inclusive_jets.size(); i++) {
00085     printf("%5u %15.8f %15.8f %15.8f\n",
00086            i, inclusive_jets[i].rap(), inclusive_jets[i].phi(),
00087            inclusive_jets[i].perp());
00088   }
00089 
00090   // simple test to avoid that the example below crashes:
00091   // make sure there is at least 2 jets above our 5 GeV
00092   if (inclusive_jets.size()<2){
00093     cout << "Please provide an event with at least 2 jets above 5 GeV" << endl;
00094     return 1;
00095   }
00096 
00097   // the sample PseudoJet that we will filter
00098   //  - the hardest jet of the event
00099   //  - the composition of the second and third hardest jets 
00100   ///   (this shows that the Filter can also be applied to a composite jet)
00101   //----------------------------------------------------------
00102   vector<PseudoJet> candidates;
00103   candidates.push_back(inclusive_jets[0]);
00104   candidates.push_back(join(inclusive_jets[1],inclusive_jets[2]));
00105 
00106 
00107   // create 5 filters
00108   //----------------------------------------------------------
00109   vector<Filter> filters;
00110   
00111   // 1.
00112   // the Cambridge/Aachen filter with Rfilt=0.3 (simpliefied version of arXiv:0802.2470)
00113   filters.push_back(Filter(JetDefinition(cambridge_algorithm, 0.3), SelectorNHardest(3)));
00114 
00115   // 2.
00116   // the Cambridge/Aachen filter with Rfilt=min(0.3, 0.5*Rbb) as in arXiv:0802.2470
00117   SharedPtr<DynamicRfilt> dynamic_Rfilt(new DynamicRfilt(0.3, 0.5));
00118   filters.push_back(Filter(dynamic_Rfilt.get(), SelectorNHardest(3)));
00119 
00120   // 3.
00121   // Filtering with a pt cut as for trimming (arXiv:0912.1342)
00122   filters.push_back(Filter(JetDefinition(kt_algorithm, 0.2), SelectorPtFractionMin(0.03)));
00123 
00124   // 4.
00125   // First example of filtering with subtraction of the background: provide rho
00126   // First, estimate the background for the given event
00127   GridMedianBackgroundEstimator bkgd(4.5, 0.55); // uses particles up to |y|=4.5
00128   bkgd.set_particles(input_particles);
00129   double rho = bkgd.rho();
00130   // Then, define the filter
00131   filters.push_back(Filter(JetDefinition(cambridge_algorithm, 0.3), SelectorNHardest(3), rho));
00132 
00133   // 5.
00134   // Second example of filtering with subtraction of the background: set a subtractor
00135   // First, define a subtractor from a background estimator
00136   Subtractor subtractor(&bkgd);
00137   // Then, define the filter
00138   Filter filt(JetDefinition(cambridge_algorithm, 0.3), SelectorNHardest(3));
00139   // Finally, tell the filter about the subtractor
00140   filt.set_subtractor(&subtractor);
00141   filters.push_back(filt);
00142 
00143 
00144   // apply the various filters to the test PseudoJet
00145   // and show the result
00146   //----------------------------------------------------------
00147 
00148   // print out original jet candidates
00149   cout << "\nOriginal jets that will be filtered: " << endl;
00150   for (vector<PseudoJet>::iterator jit=candidates.begin(); jit!=candidates.end(); jit++){
00151     const PseudoJet & c = *jit;
00152     cout << "  rap = " << c.rap() << ", phi = " << c.phi() << ", pt = " << c.perp() 
00153          << "  [" << c.description() << "]" <<  endl;
00154   }
00155 
00156   // loop on filters
00157   for (vector<Filter>::iterator it=filters.begin(); it!=filters.end(); it++){
00158     const Filter & f = *it;
00159     cout << "\nUsing filter: " << f.description() << endl;
00160     
00161     // loop on jet candidates
00162     for (vector<PseudoJet>::iterator jit=candidates.begin(); jit!=candidates.end(); jit++){
00163       const PseudoJet & c = *jit;
00164       
00165       // apply filter j to jet c      
00166       PseudoJet j = f(c);
00167       
00168       // access properties specific to the Filter
00169       //
00170       // We first make sure that the jet indeed has a structure
00171       // compatible with the result of a Filter (using
00172       // has_structure_of()), and then retrieve the pieces rejected by the
00173       // filter (using structure_of())
00174       assert(j.has_structure_of<Filter>());
00175       const Filter::StructureType & fj_struct = j.structure_of<Filter>();
00176       
00177       // write out result
00178       cout << "  rap = " << j.rap() << ", phi = " << j.phi() << ", pt = " << j.perp() 
00179            << "  [kept: " << j.pieces().size() << ", rejected: "
00180            << fj_struct.rejected().size() << "]" << endl;
00181     }
00182   }
00183 
00184   return 0;
00185 }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends