FastJet  3.1.1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Modules Pages
ClusterSequenceVoronoiArea.hh
1 //FJSTARTHEADER
2 // $Id: ClusterSequenceVoronoiArea.hh 3433 2014-07-23 08:17:03Z salam $
3 //
4 // Copyright (c) 2005-2014, 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. They are described in the original FastJet paper,
16 // hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
17 // FastJet as part of work towards a scientific publication, please
18 // quote the version you use and include a citation to the manual and
19 // optionally also to hep-ph/0512210.
20 //
21 // FastJet is distributed in the hope that it will be useful,
22 // but WITHOUT ANY WARRANTY; without even the implied warranty of
23 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 // GNU General Public License for more details.
25 //
26 // You should have received a copy of the GNU General Public License
27 // along with FastJet. If not, see <http://www.gnu.org/licenses/>.
28 //----------------------------------------------------------------------
29 //FJENDHEADER
30 
31 #ifndef __FASTJET_CLUSTERSEQUENCEVORONOIAREA_HH__
32 #define __FASTJET_CLUSTERSEQUENCEVORONOIAREA_HH__
33 
34 #include "fastjet/PseudoJet.hh"
35 #include "fastjet/AreaDefinition.hh"
36 #include "fastjet/ClusterSequenceAreaBase.hh"
37 #include <memory>
38 #include <vector>
39 
40 FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
41 
42 /// @ingroup sec_area_classes
43 /// \class ClusterSequenceVoronoiArea
44 /// Like ClusterSequence with computation of the Voronoi jet area
45 ///
46 /// Handle the computation of Voronoi jet area.
47 ///
48 /// This class should not be used directly. Rather use
49 /// ClusterSequenceArea with the appropriate AreaDefinition
51 public:
52  /// template ctor
53  /// \param pseudojet list of jets (template type)
54  /// \param jet_def jet definition
55  /// \param effective_Rfact effective radius
56  /// \param writeout_combinations ??????
57  template<class L> ClusterSequenceVoronoiArea
58  (const std::vector<L> & pseudojets,
59  const JetDefinition & jet_def,
60  const VoronoiAreaSpec & spec = VoronoiAreaSpec(),
61  const bool & writeout_combinations = false);
62 
63  /// default dtor
65 
66  /// return the area associated with the given jet
67  virtual inline double area(const PseudoJet & jet) const {
68  return _voronoi_area[jet.cluster_hist_index()];}
69 
70  /// return a 4-vector area associated with the given jet -- strictly
71  /// this is not the exact 4-vector area, but rather an approximation
72  /// made of sums of centres of all Voronoi cells in jet, each
73  /// contributing with a normalisation equal to the area of the cell
74  virtual inline PseudoJet area_4vector(const PseudoJet & jet) const {
75  return _voronoi_area_4vector[jet.cluster_hist_index()];}
76 
77  /// return the error of the area associated with the given jet
78  /// (0 by definition for a voronoi area)
79  virtual inline double area_error(const PseudoJet & /*jet*/) const {
80  return 0.0;}
81 
82  /// passive area calculator -- to be defined in the .cc file (it will do
83  /// the true hard work)
84  class VoronoiAreaCalc;
85 
86 
87 private:
88  /// initialisation of the Voronoi Area
89  void _initializeVA();
90 
91  std::vector<double> _voronoi_area; ///< vector containing the result
92  std::vector<PseudoJet> _voronoi_area_4vector; ///< vector containing approx 4-vect areas
93  VoronoiAreaCalc *_pa_calc; ///< area calculator
94  double _effective_Rfact; ///< effective radius
95 };
96 
97 
98 
99 
100 /// template constructor need to be specified in the header!
101 //----------------------------------------------------------------------
102 template<class L> ClusterSequenceVoronoiArea::ClusterSequenceVoronoiArea
103 (const std::vector<L> &pseudojets,
104  const JetDefinition &jet_def_in,
105  const VoronoiAreaSpec & spec,
106  const bool & writeout_combinations) :
107  _effective_Rfact(spec.effective_Rfact()) {
108 
109  // transfer the initial jets (type L) into our own array
110  _transfer_input_jets(pseudojets);
111 
112  // run the clustering
113  _initialise_and_run(jet_def_in,writeout_combinations);
114 
115  // the jet clustering's already been done, now worry about areas...
116  _initializeVA();
117 }
118 
119 FASTJET_END_NAMESPACE
120 
121 #endif // __FASTJET_CLUSTERSEQUENCEVORONOIAREA_HH__
virtual PseudoJet area_4vector(const PseudoJet &jet) const
return a 4-vector area associated with the given jet – strictly this is not the exact 4-vector area...
int cluster_hist_index() const
return the cluster_hist_index, intended to be used by clustering routines.
Definition: PseudoJet.hh:770
double effective_Rfact() const
return the value of effective_Rfact
virtual double area_error(const PseudoJet &) const
return the error of the area associated with the given jet (0 by definition for a voronoi area) ...
Specification for the computation of the Voronoi jet area.
base class that sets interface for extensions of ClusterSequence that provide information about the a...
virtual double area(const PseudoJet &jet) const
return the area associated with the given jet
Like ClusterSequence with computation of the Voronoi jet area.
Class to contain pseudojets, including minimal information of use to jet-clustering routines...
Definition: PseudoJet.hh:67
class that is intended to hold a full definition of the jet clusterer