Fork me on GitHub

source: git/external/fastjet/ClusterSequenceVoronoiArea.hh@ 41376ac

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 41376ac was d7d2da3, checked in by pavel <pavel@…>, 12 years ago

move branches/ModularDelphes to trunk

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