Fork me on GitHub

source: git/external/fastjet/ClusterSequenceArea.hh@ 4f728cd

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 4f728cd was 35cdc46, checked in by Pavel Demin <demin@…>, 10 years ago

upgrade FastJet to version 3.1.0-beta.1, upgrade Nsubjettiness to version 2.1.0, add SoftKiller version 1.0.0

  • Property mode set to 100644
File size: 10.7 KB
Line 
1//FJSTARTHEADER
2// $Id: ClusterSequenceArea.hh 3484 2014-07-29 21:39:39Z soyez $
3//
4// Copyright (c) 2006-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_CLUSTERSEQUENCEAREA_HH__
32#define __FASTJET_CLUSTERSEQUENCEAREA_HH__
33
34#include "fastjet/ClusterSequenceAreaBase.hh"
35#include "fastjet/ClusterSequenceActiveArea.hh"
36#include "fastjet/ClusterSequenceActiveAreaExplicitGhosts.hh"
37#include "fastjet/ClusterSequencePassiveArea.hh"
38#include "fastjet/ClusterSequenceVoronoiArea.hh"
39#include "fastjet/AreaDefinition.hh"
40
41FASTJET_BEGIN_NAMESPACE
42
43/// @ingroup area_classes
44/// \class ClusterSequenceArea
45/// General class for user to obtain ClusterSequence with additional
46/// area information.
47///
48/// Based on the area_def, it automatically dispatches the work to the
49/// appropriate actual ClusterSequenceAreaBase-derived-class to do the
50/// real work.
51class ClusterSequenceArea : public ClusterSequenceAreaBase {
52public:
53 /// main constructor
54 template<class L> ClusterSequenceArea
55 (const std::vector<L> & pseudojets,
56 const JetDefinition & jet_def_in,
57 const AreaDefinition & area_def_in) : _area_def(area_def_in) {
58 initialize_and_run_cswa(pseudojets, jet_def_in);
59 }
60
61 /// constructor with a GhostedAreaSpec
62 template<class L> ClusterSequenceArea
63 (const std::vector<L> & pseudojets,
64 const JetDefinition & jet_def_in,
65 const GhostedAreaSpec & ghost_spec) : _area_def(ghost_spec){
66 initialize_and_run_cswa(pseudojets, jet_def_in);
67 }
68
69 /// constructor with a VoronoiAreaSpec
70 template<class L> ClusterSequenceArea
71 (const std::vector<L> & pseudojets,
72 const JetDefinition & jet_def_in,
73 const VoronoiAreaSpec & voronoi_spec) : _area_def(voronoi_spec){
74 initialize_and_run_cswa(pseudojets, jet_def_in);
75 }
76
77 /// return a reference to the area definition
78 const AreaDefinition & area_def() const {return _area_def;}
79
80
81 /// return the area associated with the given jet
82 virtual double area (const PseudoJet & jet) const {
83 return _area_base->area(jet);}
84
85 /// return the error (uncertainty) associated with the determination
86 /// of the area of this jet
87 virtual double area_error (const PseudoJet & jet) const {
88 return _area_base->area_error(jet);}
89
90 /// return the 4-vector area
91 virtual PseudoJet area_4vector(const PseudoJet & jet) const {
92 return _area_base->area_4vector(jet);}
93
94 // /// return the total area, up to |y|<maxrap, that is free of jets
95 // virtual double empty_area(double maxrap) const {
96 // return _area_base->empty_area(maxrap);}
97 //
98 // /// return something similar to the number of pure ghost jets
99 // /// in the given rapidity range in an active area case.
100 // /// For the local implementation we return empty_area/(0.55 pi R^2),
101 // /// based on measured properties of ghost jets with kt and cam. Note
102 // /// that the number returned is a double.
103 // virtual double n_empty_jets(double maxrap) const {
104 // return _area_base->n_empty_jets(maxrap);
105
106 /// return the total area, corresponding to the given selector, that
107 /// is free of jets
108 ///
109 /// The selector needs to have a finite area and be applicable jet by
110 /// jet (see the BackgroundEstimator and Subtractor tools for more
111 /// advanced usage)
112 virtual double empty_area(const Selector & selector) const {
113 return _area_base->empty_area(selector);}
114
115 /// return something similar to the number of pure ghost jets
116 /// in the given rap-phi range in an active area case.
117 /// For the local implementation we return empty_area/(0.55 pi R^2),
118 /// based on measured properties of ghost jets with kt and cam. Note
119 /// that the number returned is a double.
120 ///
121 /// The selector needs to have a finite area and be applicable jet by
122 /// jet (see the BackgroundEstimator and Subtractor tools for more
123 /// advanced usage)
124 virtual double n_empty_jets(const Selector & selector) const {
125 return _area_base->n_empty_jets(selector);
126 }
127
128 /// true if a jet is made exclusively of ghosts
129 virtual bool is_pure_ghost(const PseudoJet & jet) const {
130 return _area_base->is_pure_ghost(jet);
131 }
132
133 /// true if this ClusterSequence has explicit ghosts
134 virtual bool has_explicit_ghosts() const {
135 return _area_base->has_explicit_ghosts();
136 }
137
138
139 /// overload version of what's in the ClusterSequenceAreaBase class, which
140 /// additionally checks compatibility between "selector" and region in which
141 /// ghosts are thrown.
142 ///
143 /// The selector needs to have a finite area and be applicable jet by
144 /// jet (see the BackgroundEstimator and Subtractor tools for more
145 /// advanced usage)
146 virtual void get_median_rho_and_sigma(const std::vector<PseudoJet> & all_jets,
147 const Selector & selector,
148 bool use_area_4vector,
149 double & median, double & sigma,
150 double & mean_area,
151 bool all_are_incl = false) const {
152 _warn_if_range_unsuitable(selector);
153 ClusterSequenceAreaBase::get_median_rho_and_sigma(
154 all_jets, selector, use_area_4vector,
155 median, sigma, mean_area, all_are_incl);
156 }
157
158 /// overload version of what's in the ClusterSequenceAreaBase class,
159 /// which actually just does the same thing as the base version (but
160 /// since we've overridden the 5-argument version above, we have to
161 /// override the 4-argument version too.
162 virtual void get_median_rho_and_sigma(const Selector & selector,
163 bool use_area_4vector,
164 double & median, double & sigma) const {
165 ClusterSequenceAreaBase::get_median_rho_and_sigma(selector,use_area_4vector,
166 median,sigma);
167 }
168
169 /// overload version of what's in the ClusterSequenceAreaBase class,
170 /// which actually just does the same thing as the base version (but
171 /// since we've overridden the multi-argument version above, we have to
172 /// override the 5-argument version too.
173 virtual void get_median_rho_and_sigma(const Selector & selector,
174 bool use_area_4vector,
175 double & median, double & sigma,
176 double & mean_area) const {
177 ClusterSequenceAreaBase::get_median_rho_and_sigma(selector,use_area_4vector,
178 median,sigma, mean_area);
179 }
180
181
182 /// overload version of what's in the ClusterSequenceAreaBase class, which
183 /// additionally checks compatibility between "range" and region in which
184 /// ghosts are thrown.
185 virtual void parabolic_pt_per_unit_area(double & a, double & b,
186 const Selector & selector,
187 double exclude_above=-1.0,
188 bool use_area_4vector=false) const {
189 _warn_if_range_unsuitable(selector);
190 ClusterSequenceAreaBase::parabolic_pt_per_unit_area(
191 a,b,selector, exclude_above, use_area_4vector);
192 }
193
194
195private:
196
197 /// print a warning if the range is unsuitable for the current
198 /// calculation of the area (e.g. because ghosts do not extend
199 /// far enough).
200 void _warn_if_range_unsuitable(const Selector & selector) const;
201
202 template<class L> void initialize_and_run_cswa (
203 const std::vector<L> & pseudojets,
204 const JetDefinition & jet_def);
205
206 std::auto_ptr<ClusterSequenceAreaBase> _area_base;
207 AreaDefinition _area_def;
208 static LimitedWarning _range_warnings;
209 static LimitedWarning _explicit_ghosts_repeats_warnings;
210
211};
212
213//----------------------------------------------------------------------
214template<class L> void ClusterSequenceArea::initialize_and_run_cswa(
215 const std::vector<L> & pseudojets,
216 const JetDefinition & jet_def_in)
217 {
218
219 ClusterSequenceAreaBase * _area_base_ptr;
220 switch(_area_def.area_type()) {
221 case active_area:
222 _area_base_ptr = new ClusterSequenceActiveArea(pseudojets,
223 jet_def_in,
224 _area_def.ghost_spec());
225 break;
226 case active_area_explicit_ghosts:
227 if (_area_def.ghost_spec().repeat() != 1)
228 _explicit_ghosts_repeats_warnings.warn("Requested active area with explicit ghosts with repeat != 1; only 1 set of ghosts will be used");
229 _area_base_ptr = new ClusterSequenceActiveAreaExplicitGhosts(pseudojets,
230 jet_def_in,
231 _area_def.ghost_spec());
232 break;
233 case voronoi_area:
234 _area_base_ptr = new ClusterSequenceVoronoiArea(pseudojets,
235 jet_def_in,
236 _area_def.voronoi_spec());
237 break;
238 case one_ghost_passive_area:
239 _area_base_ptr = new ClusterSequence1GhostPassiveArea(pseudojets,
240 jet_def_in,
241 _area_def.ghost_spec());
242 break;
243 case passive_area:
244 _area_base_ptr = new ClusterSequencePassiveArea(pseudojets,
245 jet_def_in,
246 _area_def.ghost_spec());
247 break;
248 default:
249 std::ostringstream err;
250 err << "Error: unrecognized area_type in ClusterSequenceArea:"
251 << _area_def.area_type();
252 throw Error(err.str());
253 //exit(-1);
254 }
255 // now copy across the information from the area base class
256 _area_base = std::auto_ptr<ClusterSequenceAreaBase>(_area_base_ptr);
257 transfer_from_sequence(*_area_base);
258}
259
260FASTJET_END_NAMESPACE
261
262#endif // __FASTJET_CLUSTERSEQUENCEAREA_HH__
263
264
Note: See TracBrowser for help on using the repository browser.