31 #include "fastjet/ClusterSequenceActiveAreaExplicitGhosts.hh" 36 FASTJET_BEGIN_NAMESPACE
39 typedef ClusterSequenceActiveAreaExplicitGhosts ClustSeqActAreaEG;
42 LimitedWarning ClustSeqActAreaEG::_warnings;
46 void ClustSeqActAreaEG::_add_ghosts (
47 const GhostedAreaSpec & ghost_spec) {
50 ghost_spec.add_ghosts(_jets);
53 for (
unsigned i = _initial_hard_n; i < _jets.size(); i++) {
55 _is_pure_ghost.push_back(
true);
59 _ghost_area = ghost_spec.actual_ghost_area();
60 _n_ghosts = ghost_spec.n_ghosts();
66 double ClustSeqActAreaEG::area (
const PseudoJet & jet)
const {
73 double ClustSeqActAreaEG::total_area ()
const {
74 return _n_ghosts * _ghost_area;
85 bool ClustSeqActAreaEG::is_pure_ghost(
const PseudoJet & jet)
const 91 bool ClustSeqActAreaEG::is_pure_ghost(
int hist_ix)
const 93 return hist_ix >= 0 ? _is_pure_ghost[hist_ix] :
false;
97 double ClustSeqActAreaEG::empty_area(
const Selector & selector)
const {
100 throw Error(
"ClusterSequenceActiveAreaExplicitGhosts: empty area can only be computed from selectors applying jet by jet");
103 vector<PseudoJet> unclust = unclustered_particles();
104 double area_local = 0.0;
105 for (
unsigned iu = 0; iu < unclust.size(); iu++) {
106 if (is_pure_ghost(unclust[iu]) && selector.
pass(unclust[iu])) {
107 area_local += _ghost_area;
115 void ClustSeqActAreaEG::_post_process() {
119 _max_ghost_perp2 = 0.0;
120 for (
int i = 0; i < _initial_n; i++) {
121 if (_is_pure_ghost[i] && _jets[i].perp2() > _max_ghost_perp2)
122 _max_ghost_perp2 = _jets[i].perp2();
126 double danger_ratio = numeric_limits<double>::epsilon();
127 danger_ratio = danger_ratio * danger_ratio;
128 _has_dangerous_particles =
false;
129 for (
int i = 0; i < _initial_n; i++) {
130 if (!_is_pure_ghost[i] &&
131 danger_ratio * _jets[i].perp2() <= _max_ghost_perp2) {
132 _has_dangerous_particles =
true;
137 if (_has_dangerous_particles) _warnings.warn(
"ClusterSequenceActiveAreaExplicitGhosts: \n ghosts not sufficiently soft wrt some of the input particles\n a common cause is (unphysical?) input particles with pt=0 but finite rapidity");
140 _areas.resize(_history.size());
141 _area_4vectors.resize(_history.size());
142 _is_pure_ghost.resize(_history.size());
160 for (
int i = 0; i < _initial_n; i++) {
161 if (_is_pure_ghost[i]) {
162 _areas[i] = _ghost_area;
171 _area_4vectors[i].reset_momentum(_jets[i]);
172 _area_4vectors[i] *= (_ghost_area/_jets[i].perp());
175 _area_4vectors[i] =
PseudoJet(0.0,0.0,0.0,0.0);
182 for (
unsigned i = _initial_n; i < _history.size(); i++) {
183 if (_history[i].parent2 == BeamJet) {
184 _is_pure_ghost[i] = _is_pure_ghost[_history[i].parent1];
185 _areas[i] = _areas[_history[i].parent1];
186 _area_4vectors[i] = _area_4vectors[_history[i].parent1];
188 _is_pure_ghost[i] = _is_pure_ghost[_history[i].parent1] &&
189 _is_pure_ghost[_history[i].parent2] ;
190 _areas[i] = _areas[_history[i].parent1] +
191 _areas[_history[i].parent2] ;
192 _jet_def.recombiner()->recombine(_area_4vectors[_history[i].parent1],
193 _area_4vectors[_history[i].parent2],
203 FASTJET_END_NAMESPACE
bool applies_jet_by_jet() const
returns true if this can be applied jet by jet
bool pass(const PseudoJet &jet) const
return true if the jet passes the selection
base class corresponding to errors that can be thrown by FastJet
int cluster_hist_index() const
return the cluster_hist_index, intended to be used by clustering routines.
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
Class to contain pseudojets, including minimal information of use to jet-clustering routines...