Fork me on GitHub

source: git/external/fastjet/tools/BackgroundEstimatorBase.hh@ a0c065d

ImprovedOutputFile Timing dual_readout llp
Last change on this file since a0c065d 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: 9.3 KB
Line 
1#ifndef __FASTJET_BACKGROUND_ESTIMATOR_BASE_HH__
2#define __FASTJET_BACKGROUND_ESTIMATOR_BASE_HH__
3
4//FJSTARTHEADER
5// $Id: BackgroundEstimatorBase.hh 3516 2014-08-01 14:07:58Z salam $
6//
7// Copyright (c) 2005-2014, Matteo Cacciari, Gavin P. Salam and Gregory Soyez
8//
9//----------------------------------------------------------------------
10// This file is part of FastJet.
11//
12// FastJet is free software; you can redistribute it and/or modify
13// it under the terms of the GNU General Public License as published by
14// the Free Software Foundation; either version 2 of the License, or
15// (at your option) any later version.
16//
17// The algorithms that underlie FastJet have required considerable
18// development. They are described in the original FastJet paper,
19// hep-ph/0512210 and in the manual, arXiv:1111.6097. If you use
20// FastJet as part of work towards a scientific publication, please
21// quote the version you use and include a citation to the manual and
22// optionally also to hep-ph/0512210.
23//
24// FastJet is distributed in the hope that it will be useful,
25// but WITHOUT ANY WARRANTY; without even the implied warranty of
26// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
27// GNU General Public License for more details.
28//
29// You should have received a copy of the GNU General Public License
30// along with FastJet. If not, see <http://www.gnu.org/licenses/>.
31//----------------------------------------------------------------------
32//FJENDHEADER
33
34#include <fastjet/ClusterSequenceAreaBase.hh>
35#include <fastjet/FunctionOfPseudoJet.hh>
36#include <fastjet/Selector.hh>
37#include <fastjet/Error.hh>
38#include <iostream>
39
40FASTJET_BEGIN_NAMESPACE // defined in fastjet/internal/base.hh
41
42
43/// @ingroup tools_background
44/// \class BackgroundEstimatorBase
45///
46/// Abstract base class that provides the basic interface for classes
47/// that estimate levels of background radiation in hadron and
48/// heavy-ion collider events.
49///
50class BackgroundEstimatorBase {
51public:
52 /// @name constructors and destructors
53 //\{
54 //----------------------------------------------------------------
55 BackgroundEstimatorBase() : _rescaling_class(0) {}
56 //\}
57
58 /// a default virtual destructor that does nothing
59 virtual ~BackgroundEstimatorBase() {}
60
61
62 /// @name setting a new event
63 //\{
64 //----------------------------------------------------------------
65
66 /// tell the background estimator that it has a new event, composed
67 /// of the specified particles.
68 virtual void set_particles(const std::vector<PseudoJet> & particles) = 0;
69
70 //\}
71
72 /// @name retrieving fundamental information
73 //\{
74 //----------------------------------------------------------------
75
76 /// get rho, the background density per unit area
77 virtual double rho() const = 0;
78
79 /// get sigma, the background fluctuations per unit area; must be
80 /// multipled by sqrt(area) to get fluctuations for a region of a
81 /// given area.
82 virtual double sigma() const {
83 throw Error("sigma() not supported for this Background Estimator");
84 }
85
86 /// get rho, the background density per unit area, locally at the
87 /// position of a given jet. Note that this is not const, because a
88 /// user may then wish to query other aspects of the background that
89 /// could depend on the position of the jet last used for a rho(jet)
90 /// determination.
91 virtual double rho(const PseudoJet & jet) = 0;
92
93 /// get sigma, the background fluctuations per unit area, locally at
94 /// the position of a given jet. As for rho(jet), it is non-const.
95 virtual double sigma(const PseudoJet & /*jet*/) {
96 throw Error("sigma(jet) not supported for this Background Estimator");
97 }
98
99 /// returns true if this background estimator has support for
100 /// determination of sigma
101 virtual bool has_sigma() {return false;}
102
103 //----------------------------------------------------------------
104 // now do the same thing for rho_m and sigma_m
105
106 /// returns rho_m, the purely longitudinal, particle-mass-induced
107 /// component of the background density per unit area
108 virtual double rho_m() const{
109 throw Error("rho_m() not supported for this Background Estimator");
110 }
111
112 /// returns sigma_m, a measure of the fluctuations in the purely
113 /// longitudinal, particle-mass-induced component of the background
114 /// density per unit area; must be multipled by sqrt(area) to get
115 /// fluctuations for a region of a given area.
116 virtual double sigma_m() const {
117 throw Error("sigma_m() not supported for this Background Estimator");
118 }
119
120 /// Returns rho_m locally at the jet position. As for rho(jet), it is non-const.
121 virtual double rho_m(const PseudoJet & /*jet*/){
122 throw Error("rho_m(jet) not supported for this Background Estimator");
123 }
124
125 /// Returns sigma_m locally at the jet position. As for rho(jet), it is non-const.
126 virtual double sigma_m(const PseudoJet & /*jet*/) {
127 throw Error("sigma_m(jet) not supported for this Background Estimator");
128 }
129
130 /// Returns true if this background estimator has support for
131 /// determination of rho_m.
132 ///
133 /// Note that support for sigma_m is automatic is one has sigma and
134 /// rho_m support.
135 virtual bool has_rho_m() const {return false;}
136 //\}
137
138
139 /// @name configuring the behaviour
140 //\{
141 //----------------------------------------------------------------
142
143 /// Set a pointer to a class that calculates the rescaling factor as
144 /// a function of the jet (position). Note that the rescaling factor
145 /// is used both in the determination of the "global" rho (the pt/A
146 /// of each jet is divided by this factor) and when asking for a
147 /// local rho (the result is multiplied by this factor).
148 ///
149 /// The BackgroundRescalingYPolynomial class can be used to get a
150 /// rescaling that depends just on rapidity.
151 ///
152 /// There is currently no support for different rescaling classes
153 /// for rho and rho_m determinations.
154 virtual void set_rescaling_class(const FunctionOfPseudoJet<double> * rescaling_class_in) { _rescaling_class = rescaling_class_in; }
155
156 /// return the pointer to the jet density class
157 const FunctionOfPseudoJet<double> * rescaling_class() const{
158 return _rescaling_class;
159 }
160
161 //\}
162
163 /// @name description
164 //\{
165 //----------------------------------------------------------------
166
167 /// returns a textual description of the background estimator
168 virtual std::string description() const = 0;
169
170 //\}
171
172protected:
173 /// @name helpers for derived classes
174 ///
175 /// Note that these helpers are related to median-based estimation
176 /// of the background, so there is no guarantee that they will
177 /// remain in this base class in the long term
178 //\{
179 //----------------------------------------------------------------
180
181 /// given a quantity in a vector (e.g. pt_over_area) and knowledge
182 /// about the number of empty jets, calculate the median and
183 /// stand_dev_if_gaussian (roughly from the 16th percentile)
184 ///
185 /// If do_fj2_calculation is set to true then this performs FastJet
186 /// 2.X estimation of the standard deviation, which has a spurious
187 /// offset in the limit of a small number of jets.
188 void _median_and_stddev(const std::vector<double> & quantity_vector,
189 double n_empty_jets,
190 double & median,
191 double & stand_dev_if_gaussian,
192 bool do_fj2_calculation = false
193 ) const;
194
195 /// computes a percentile of a given _sorted_ vector
196 /// \param sorted_quantity_vector the vector contains the data sample
197 /// \param percentile the percentile (defined between 0 and 1) to compute
198 /// \param nempty an additional number of 0's
199 /// (considered at the beginning of
200 /// the quantity vector)
201 /// \param do_fj2_calculation carry out the calculation as it
202 /// was done in fj2 (suffers from "edge effects")
203 double _percentile(const std::vector<double> & sorted_quantity_vector,
204 const double percentile,
205 const double nempty=0.0,
206 const bool do_fj2_calculation = false) const;
207
208 //\}
209
210 const FunctionOfPseudoJet<double> * _rescaling_class;
211 static LimitedWarning _warnings_empty_area;
212};
213
214
215
216//----------------------------------------------------------------------
217/// @ingroup tools_background
218/// A background rescaling that is a simple polynomial in y
219class BackgroundRescalingYPolynomial : public FunctionOfPseudoJet<double> {
220public:
221 /// construct a background rescaling polynomial of the form
222 /// a0 + a1*y + a2*y^2 + a3*y^3 + a4*y^4
223 ///
224 /// The following values give a reasonable reproduction of the
225 /// Pythia8 tune 4C background shape for pp collisions at
226 /// sqrt(s)=7TeV:
227 ///
228 /// - a0 = 1.157
229 /// - a1 = 0
230 /// - a2 = -0.0266
231 /// - a3 = 0
232 /// - a4 = 0.000048
233 ///
234 BackgroundRescalingYPolynomial(double a0=1,
235 double a1=0,
236 double a2=0,
237 double a3=0,
238 double a4=0) : _a0(a0), _a1(a1), _a2(a2), _a3(a3), _a4(a4) {}
239
240 /// return the rescaling factor associated with this jet
241 virtual double result(const PseudoJet & jet) const;
242private:
243 double _a0, _a1, _a2, _a3, _a4;
244};
245
246
247
248
249
250FASTJET_END_NAMESPACE
251
252#endif // __BACKGROUND_ESTIMATOR_BASE_HH__
253
Note: See TracBrowser for help on using the repository browser.