FastJet 3.0beta1
ClusterSequenceVoronoiArea.hh
00001 //STARTHEADER
00002 // $Id: ClusterSequenceVoronoiArea.hh 1761 2010-09-16 10:43:18Z soyez $
00003 //
00004 // Copyright (c) 2005-2007, Matteo Cacciari, Gavin Salam and Gregory Soyez
00005 //
00006 //----------------------------------------------------------------------
00007 // This file is part of FastJet.
00008 //
00009 //  FastJet is free software; you can redistribute it and/or modify
00010 //  it under the terms of the GNU General Public License as published by
00011 //  the Free Software Foundation; either version 2 of the License, or
00012 //  (at your option) any later version.
00013 //
00014 //  The algorithms that underlie FastJet have required considerable
00015 //  development and are described in hep-ph/0512210. If you use
00016 //  FastJet as part of work towards a scientific publication, please
00017 //  include a citation to the FastJet paper.
00018 //
00019 //  FastJet is distributed in the hope that it will be useful,
00020 //  but WITHOUT ANY WARRANTY; without even the implied warranty of
00021 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00022 //  GNU General Public License for more details.
00023 //
00024 //  You should have received a copy of the GNU General Public License
00025 //  along with FastJet; if not, write to the Free Software
00026 //  Foundation, Inc.:
00027 //      59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00028 //----------------------------------------------------------------------
00029 //ENDHEADER
00030 
00031 #ifndef __FASTJET_CLUSTERSEQUENCEVORONOIAREA_HH__
00032 #define __FASTJET_CLUSTERSEQUENCEVORONOIAREA_HH__
00033 
00034 #include "fastjet/PseudoJet.hh"
00035 #include "fastjet/AreaDefinition.hh"
00036 #include "fastjet/ClusterSequenceAreaBase.hh"
00037 #include <memory>
00038 #include <vector>
00039 
00040 FASTJET_BEGIN_NAMESPACE      // defined in fastjet/internal/base.hh
00041 
00042 /// @ingroup sec_area_classes
00043 /// \class ClusterSequenceVoronoiArea
00044 /// Like ClusterSequence with computation of the Voronoi jet area
00045 ///
00046 /// Handle the computation of Voronoi jet area.
00047 ///
00048 /// This class should not be used directly. Rather use
00049 /// ClusterSequenceArea with the appropriate AreaDefinition
00050 class ClusterSequenceVoronoiArea : public ClusterSequenceAreaBase {
00051 public:
00052   /// template ctor
00053   /// \param pseudojet              list of jets (template type)
00054   /// \param jet_def                jet definition
00055   /// \param effective_Rfact        effective radius
00056   /// \param writeout_combinations  ??????
00057   template<class L> ClusterSequenceVoronoiArea
00058          (const std::vector<L> & pseudojets, 
00059           const JetDefinition & jet_def,
00060           const VoronoiAreaSpec & spec = VoronoiAreaSpec(),
00061           const bool & writeout_combinations = false);
00062   
00063   /// default dtor
00064   ~ClusterSequenceVoronoiArea();
00065 
00066   /// return the area associated with the given jet
00067   virtual inline double area(const PseudoJet & jet) const {
00068     return _voronoi_area[jet.cluster_hist_index()];};
00069 
00070   /// return a 4-vector area associated with the given jet -- stricly
00071   /// this is not the exact 4-vector area, but rather an approximation
00072   /// made of sums of centres of all Voronoi cells in jet, each
00073   /// contributing with a normalisation equal to the area of the cell
00074   virtual inline PseudoJet area_4vector(const PseudoJet & jet) const {
00075     return _voronoi_area_4vector[jet.cluster_hist_index()];};
00076 
00077   /// return the error of the area associated with the given jet
00078   /// (0 by definition for a voronoi area)
00079   virtual inline double area_error(const PseudoJet & jet) const {
00080     return 0.0;};
00081 
00082   /// passive area calculator -- to be defined in the .cc file (it will do
00083   /// the true hard work)
00084   class VoronoiAreaCalc; 
00085   
00086 
00087 private:  
00088   /// initialisation of the Voronoi Area
00089   void _initializeVA();
00090 
00091   std::vector<double> _voronoi_area;  ///< vector containing the result
00092   std::vector<PseudoJet> _voronoi_area_4vector; ///< vector containing approx 4-vect areas
00093   VoronoiAreaCalc *_pa_calc;          ///< area calculator
00094   double _effective_Rfact;            ///< effective radius
00095 };
00096 
00097 
00098 
00099 
00100 /// template constructor need to be specified in the header!
00101 //----------------------------------------------------------------------
00102 template<class L> ClusterSequenceVoronoiArea::ClusterSequenceVoronoiArea
00103 (const std::vector<L> &pseudojets, 
00104  const JetDefinition &jet_def,
00105  const VoronoiAreaSpec & spec,
00106  const bool & writeout_combinations) :
00107   _effective_Rfact(spec.effective_Rfact()) {
00108 
00109   // transfer the initial jets (type L) into our own array
00110   _transfer_input_jets(pseudojets);
00111 
00112   // run the clustering
00113   _initialise_and_run(jet_def,writeout_combinations);
00114 
00115   // the jet clustering's already been done, now worry about areas...
00116   _initializeVA();
00117 }
00118 
00119 FASTJET_END_NAMESPACE
00120 
00121 #endif // __FASTJET_CLUSTERSEQUENCEVORONOIAREA_HH__
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends