Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:57:18

0001 #ifndef __FASTJET_NNBASE_HH__
0002 #define __FASTJET_NNBASE_HH__
0003 
0004 //FJSTARTHEADER
0005 // $Id$
0006 //
0007 // Copyright (c) 2016-2021, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
0008 //
0009 //----------------------------------------------------------------------
0010 // This file is part of FastJet.
0011 //
0012 //  FastJet is free software; you can redistribute it and/or modify
0013 //  it under the terms of the GNU General Public License as published by
0014 //  the Free Software Foundation; either version 2 of the License, or
0015 //  (at your option) any later version.
0016 //
0017 //  The algorithms that underlie FastJet have required considerable
0018 //  development. They are described in the original FastJet paper,
0019 //  hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
0020 //  FastJet as part of work towards a scientific publication, please
0021 //  quote the version you use and include a citation to the manual and
0022 //  optionally also to hep-ph/0512210.
0023 //
0024 //  FastJet is distributed in the hope that it will be useful,
0025 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
0026 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
0027 //  GNU General Public License for more details.
0028 //
0029 //  You should have received a copy of the GNU General Public License
0030 //  along with FastJet. If not, see <http://www.gnu.org/licenses/>.
0031 //----------------------------------------------------------------------
0032 //FJENDHEADER
0033 
0034 #include<fastjet/ClusterSequence.hh>
0035 
0036 
0037 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
0038 
0039 /// @ingroup advanced_usage
0040 /// \class _NoInfo
0041 /// internal dummy class, used as a default template argument
0042 class _NoInfo {};
0043 
0044 /// @ingroup advanced_usage
0045 /// \class NNInfo
0046 ///
0047 /// internal helper template class to facilitate initialisation of a
0048 /// BJ with a PseudoJet and extra information. Implementations of
0049 /// NN-based clustering do not need to explicitly use or refer to
0050 /// this class!
0051 template<class I> class NNInfo {
0052 public:
0053   NNInfo()         : _info(NULL) {}
0054   NNInfo(I * info) : _info(info) {}
0055   template<class BJ> void init_jet(BJ * briefjet, const fastjet::PseudoJet & jet, int index) { briefjet->init(jet, index, _info);}
0056 private:
0057   I * _info;
0058 };
0059 
0060 /// @ingroup advanced_usage Internal helper specialisation of NNInfo
0061 /// for cases where there is no extra info
0062 template<> class NNInfo<_NoInfo>  {
0063 public:
0064   NNInfo()           {}
0065   NNInfo(_NoInfo * ) {}
0066   template<class BJ> void init_jet(BJ * briefjet, const fastjet::PseudoJet & jet, int index) { briefjet->init(jet, index);}
0067 };
0068 
0069 
0070 //----------------------------------------------------------------------
0071 /// @ingroup advanced_usage
0072 /// \class NNBase
0073 /// Helps solve closest pair problems with generic interparticle and
0074 /// particle-beam distances.
0075 ///
0076 /// \section Description Description and derived classes:
0077 ///
0078 ///   This is an abstract base class which defines the interface for
0079 ///   several classes that help carry out nearest-neighbour
0080 ///   clustering:
0081 ///  
0082 ///    - NNH          provides an implementation for generic measures,
0083 ///  
0084 ///    - NNFJN2Plain  provides an implementation for distances
0085 ///                   satisfying the FastJet lemma i.e. distances for
0086 ///                   which the minimum dij has the property that i is
0087 ///                   the geometrical nearest neighbour of j, or vice
0088 ///                   versa. I.e. the distance can be factorised in a
0089 ///                   momentum factor and a geometric piece. This is
0090 ///                   based on the fastjet N2Plain clustering strategy
0091 ///  
0092 ///    - NNFJN2Tiled  is a tiled version of NNFJN2Plain (based on the
0093 ///                   N2Tiled FastJet clustering strategy). Like
0094 ///                   NNPlain2 it applies to distance measures that
0095 ///                   satisfy the FastJet lemma, with the additional
0096 ///                   restriction that: (a) the underlying geometry
0097 ///                   should be cylindrical (e.g. rapidity--azimuth)
0098 ///                   and (b) the search for the geometric nearest
0099 ///                   neighbour of each particle can be limited to
0100 ///                   that particle's tile and its neighbouring tiles.
0101 ///
0102 /// If you can use NNFJN2Plain it will usually be faster than
0103 /// NNH. NNFJN2Tiled, where it can be used, will be faster for
0104 /// multiplicities above a few tens of particles.
0105 ///
0106 ///   NOTE: IN ALL CASES, THE DISTANCE MUST BE SYMMETRIC (dij=dji)!!!
0107 ///
0108 /// \section BJ Underlying BriefJet (BJ) class:
0109 /// 
0110 ///   All derived classes must be templated with a BriefJet (BJ)
0111 ///   class --- BJ should basically cache the minimal amount of
0112 ///   information that is needed to efficiently calculate
0113 ///   interparticle distances and particle-beam distances.
0114 ///   
0115 ///   This class can be used with or without an extra "Information"
0116 ///   template, i.e. `NN*<BJ>` or `NN*<BJ,I>`. Accordingly BJ must provide
0117 ///   one of the two following init functions:
0118 ///
0119 ///   \code
0120 ///     void  BJ::init(const PseudoJet & jet);            // initialise with a PseudoJet
0121 ///     void  BJ::init(const PseudoJet & jet, I * info);  // initialise with a PseudoJet + info
0122 ///   \endcode
0123 ///
0124 ///   where info might be a pointer to a class that contains, e.g.,
0125 ///   information about R, or other parameters of the jet algorithm
0126 ///   
0127 ///   The BJ then provides information about interparticle and
0128 ///   particle-beam distances. The exact requirements depend on
0129 ///   whether you use NNH, NNFJN2Plain or NNFJN2Tiled. (See the
0130 ///   corresponding classes for details).
0131 ///
0132 ///
0133 /// \section Workflow Workflow:
0134 ///
0135 ///   In all cases, the usage of NNBase classes works as follows:
0136 ///
0137 ///   First, from the list of particles, create an `NN*<BJ>`
0138 ///   object of the appropriate type with the appropriate BJ class
0139 ///   (and optional extra info).
0140 ///
0141 ///   Then, cluster using a loop like this (assuming a FastJet plugin)
0142 ///
0143 ///   \code
0144 ///     while (njets > 0) {
0145 ///       int i, j, k;
0146 ///       // get the i and j that minimize the distance
0147 ///       double dij = nn.dij_min(i, j);  
0148 ///     
0149 ///       // do the appropriate recombination and update the nn
0150 ///       if (j >= 0) {    // interparticle recombination
0151 ///         cs.plugin_record_ij_recombination(i, j, dij, k);
0152 ///         nn.merge_jets(i, j, cs.jets()[k], k); 
0153 ///       } else {         // bbeam recombination
0154 ///         double diB = cs.jets()[i].E()*cs.jets()[i].E(); // get new diB
0155 ///         cs.plugin_record_iB_recombination(i, diB);
0156 ///         nn.remove_jet(i);
0157 ///       }
0158 ///       njets--;
0159 ///     }
0160 ///   \endcode
0161 ///
0162 /// For an example of how the NNH<BJ> class is used, see the
0163 /// JadePlugin or EECambridgePlugin.
0164 template<class I = _NoInfo> class NNBase : public NNInfo<I> {
0165 public:
0166   /// Default constructor
0167   NNBase() {}
0168   /// Constuctor with additional Info 
0169   NNBase(I * info) : NNInfo<I>(info) {}
0170 
0171   /// initialisation from a given list of particles
0172   virtual void start(const std::vector<PseudoJet> & jets) = 0;
0173   
0174   /// returns the dij_min and indices iA, iB, for the corresponding jets.
0175   /// If iB < 0 then iA recombines with the beam
0176   virtual double dij_min(int & iA, int & iB) = 0;
0177 
0178   /// removes the jet pointed to by index iA
0179   virtual void remove_jet(int iA) = 0;
0180 
0181   /// merges the jets pointed to by indices A and B and replaces them with
0182   /// jet, assigning it an index jet_index.
0183   virtual void merge_jets(int iA, int iB, const PseudoJet & jet, int jet_index) =  0;
0184 
0185   virtual ~NNBase() {};
0186 };
0187 
0188 
0189 FASTJET_END_NAMESPACE      // defined in fastjet/internal/base.hh
0190 
0191 
0192 #endif // __FASTJET_NNBASE_HH__