31 #include "fastjet/PseudoJet.hh" 32 #include "fastjet/ClusterSequence.hh" 33 #include "fastjet/ClusterSequenceActiveArea.hh" 34 #include "fastjet/ClusterSequenceActiveAreaExplicitGhosts.hh" 42 FASTJET_BEGIN_NAMESPACE
53 void ClusterSequenceActiveArea::_initialise_and_run_AA (
56 const bool & writeout_combinations) {
58 bool continue_running;
59 _initialise_AA(jet_def_in, ghost_spec, writeout_combinations, continue_running);
60 if (continue_running) {
62 _postprocess_AA(ghost_spec);
67 void ClusterSequenceActiveArea::_resize_and_zero_AA () {
69 _average_area.resize(2*_jets.size()); _average_area = 0.0;
70 _average_area2.resize(2*_jets.size()); _average_area2 = 0.0;
71 _average_area_4vector.resize(2*_jets.size());
72 _average_area_4vector =
PseudoJet(0.0,0.0,0.0,0.0);
73 _non_jet_area = 0.0; _non_jet_area2 = 0.0; _non_jet_number=0.0;
77 void ClusterSequenceActiveArea::_initialise_AA (
80 const bool & writeout_combinations,
81 bool & continue_running)
85 _ghost_spec_repeat = ghost_spec.repeat();
88 _resize_and_zero_AA();
91 _maxrap_for_area = ghost_spec.ghost_maxrap();
92 _safe_rap_for_area = _maxrap_for_area - jet_def_in.R();
101 if (ghost_spec.repeat() <= 0) {
102 _initialise_and_run(jet_def_in, writeout_combinations);
103 continue_running =
false;
108 _decant_options(jet_def_in, writeout_combinations);
112 _fill_initial_history();
115 _has_dangerous_particles =
false;
117 continue_running =
true;
122 void ClusterSequenceActiveArea::_run_AA (
const GhostedAreaSpec & ghost_spec) {
124 vector<PseudoJet> input_jets(_jets);
127 vector<int> unique_tree;
130 for (
int irepeat = 0; irepeat < ghost_spec.repeat(); irepeat++) {
133 jet_def(), ghost_spec);
139 _transfer_ghost_free_history(clust_seq);
141 unique_tree = unique_history_order();
145 _transfer_areas(unique_tree, clust_seq);
152 void ClusterSequenceActiveArea::_postprocess_AA (
const GhostedAreaSpec & ghost_spec) {
153 _average_area /= ghost_spec.repeat();
154 _average_area2 /= ghost_spec.repeat();
155 if (ghost_spec.repeat() > 1) {
158 const double tmp = ghost_spec.repeat()-1;
159 _average_area2 = sqrt(abs(_average_area2 - _average_area*_average_area)/tmp);
161 _average_area2 = 0.0;
164 _non_jet_area /= ghost_spec.repeat();
165 _non_jet_area2 /= ghost_spec.repeat();
166 _non_jet_area2 = sqrt(abs(_non_jet_area2 - _non_jet_area*_non_jet_area)/
167 ghost_spec.repeat());
168 _non_jet_number /= ghost_spec.repeat();
173 for (
unsigned i = 0; i < _average_area_4vector.size(); i++) {
174 _average_area_4vector[i] = (1.0/ghost_spec.repeat()) * _average_area_4vector[i];
276 double ClusterSequenceActiveArea::pt_per_unit_area(
279 vector<PseudoJet> incl_jets = inclusive_jets();
280 vector<double> pt_over_areas;
282 for (
unsigned i = 0; i < incl_jets.size(); i++) {
283 if (abs(incl_jets[i].rap()) < _safe_rap_for_area) {
285 if ( strat == median_4vector ) {
286 this_area = area_4vector(incl_jets[i]).perp();
288 this_area = area(incl_jets[i]);
290 pt_over_areas.push_back(incl_jets[i].perp()/this_area);
295 if (pt_over_areas.size() == 0) {
return 0.0;}
300 sort(pt_over_areas.begin(), pt_over_areas.end());
301 double non_ghost_median_ratio = pt_over_areas[pt_over_areas.size()/2];
306 double nj_median_pos = (pt_over_areas.size()-1 - _non_jet_number)/2.0;
307 double nj_median_ratio;
308 if (nj_median_pos >= 0 && pt_over_areas.size() > 1) {
309 int int_nj_median = int(nj_median_pos);
311 pt_over_areas[int_nj_median] * (int_nj_median+1-nj_median_pos)
312 + pt_over_areas[int_nj_median+1] * (nj_median_pos - int_nj_median);
314 nj_median_ratio = 0.0;
319 double pt_sum = 0.0, pt_sum_with_cut = 0.0;
320 double area_sum = _non_jet_area, area_sum_with_cut = _non_jet_area;
321 double ratio_sum = 0.0;
322 double ratio_n = _non_jet_number;
323 for (
unsigned i = 0; i < incl_jets.size(); i++) {
324 if (abs(incl_jets[i].rap()) < _safe_rap_for_area) {
326 if ( strat == median_4vector ) {
327 this_area = area_4vector(incl_jets[i]).perp();
329 this_area = area(incl_jets[i]);
331 pt_sum += incl_jets[i].perp();
332 area_sum += this_area;
333 double ratio = incl_jets[i].perp()/this_area;
334 if (ratio < range*nj_median_ratio) {
335 pt_sum_with_cut += incl_jets[i].perp();
336 area_sum_with_cut += this_area;
337 ratio_sum += ratio; ratio_n++;
343 double trunc_sum = 0, trunc_sumsqr = 0;
344 vector<double> means(pt_over_areas.size()), sd(pt_over_areas.size());
345 for (
unsigned i = 0; i < pt_over_areas.size() ; i++ ) {
346 double ratio = pt_over_areas[i];
348 trunc_sumsqr += ratio*ratio;
349 means[i] = trunc_sum / (i+1);
350 sd[i] = sqrt(abs(means[i]*means[i] - trunc_sumsqr/(i+1)));
351 cerr <<
"i, means, sd: " <<i<<
", "<< means[i] <<
", "<<sd[i]<<
", "<<
352 sd[i]/sqrt(i+1.0)<<endl;
354 cout <<
"-----------------------------------"<<endl;
355 for (
unsigned i = 0; i <= pt_over_areas.size()/2 ; i++ ) {
356 cout <<
"Median "<< i <<
" = " << pt_over_areas[i]<<endl;
358 cout <<
"Number of non-jets: "<<_non_jet_number<<endl;
359 cout <<
"Area of non-jets: "<<_non_jet_area<<endl;
360 cout <<
"Default median position: " << (pt_over_areas.size()-1)/2.0<<endl;
361 cout <<
"NJ median position: " << nj_median_pos <<endl;
362 cout <<
"NJ median value: " << nj_median_ratio <<endl;
369 return nj_median_ratio;
370 case non_ghost_median:
371 return non_ghost_median_ratio;
372 case pttot_over_areatot:
373 return pt_sum / area_sum;
374 case pttot_over_areatot_cut:
375 return pt_sum_with_cut / area_sum_with_cut;
377 return ratio_sum/ratio_n;
379 return nj_median_ratio;
445 double ClusterSequenceActiveArea::empty_area(
const Selector & selector)
const {
448 throw Error(
"ClusterSequenceActiveArea: empty area can only be computed from selectors applying jet by jet");
453 for (
unsigned i = 0; i < _ghost_jets.size(); i++) {
454 if (selector.
pass(_ghost_jets[i])) {
455 empty += _ghost_jets[i].area;
459 for (
unsigned i = 0; i < _unclustered_ghosts.size(); i++) {
460 if (selector.
pass(_unclustered_ghosts[i])) {
461 empty += _unclustered_ghosts[i].area;
464 empty /= _ghost_spec_repeat;
469 double ClusterSequenceActiveArea::n_empty_jets(
const Selector & selector)
const {
470 _check_selector_good_for_median(selector);
473 for (
unsigned i = 0; i < _ghost_jets.size(); i++) {
474 if (selector.
pass(_ghost_jets[i])) inrange++;
476 inrange /= _ghost_spec_repeat;
483 void ClusterSequenceActiveArea::_transfer_ghost_free_history(
486 const vector<history_element> & gs_history = ghosted_seq.
history();
487 vector<int> gs2self_hist_map(gs_history.size());
496 while (igs < gs_history.size() && gs_history[igs].parent1 == InexistentParent) {
499 gs2self_hist_map[igs] = iself++;
501 gs2self_hist_map[igs] = Invalid;
508 assert(iself == _history.size());
514 if (igs == gs_history.size())
return;
520 gs2self_hist_map[igs] = Invalid;
526 bool parent1_is_ghost = ghosted_seq.
is_pure_ghost(gs_hist_el.parent1);
532 if (parent1_is_ghost && !parent2_is_ghost && gs_hist_el.
parent2 >= 0) {
533 gs2self_hist_map[igs] = gs2self_hist_map[gs_hist_el.
parent2];
536 if (!parent1_is_ghost && parent2_is_ghost) {
537 gs2self_hist_map[igs] = gs2self_hist_map[gs_hist_el.parent1];
544 gs2self_hist_map[igs] = _history.size();
547 int jet_i = _history[gs2self_hist_map[gs_hist_el.parent1]].jetp_index;
548 int jet_j = _history[gs2self_hist_map[gs_hist_el.
parent2]].jetp_index;
550 _do_ij_recombination_step(jet_i, jet_j, gs_hist_el.
dij, newjet_k);
553 assert(gs_history[igs].parent2 == BeamJet);
555 gs2self_hist_map[igs] = _history.size();
557 _do_iB_recombination_step(
558 _history[gs2self_hist_map[gs_hist_el.parent1]].jetp_index,
561 }
while (++igs < gs_history.size());
566 void ClusterSequenceActiveArea::_transfer_areas(
567 const vector<int> & unique_hist_order,
570 const vector<history_element> & gs_history = ghosted_seq.
history();
571 const vector<PseudoJet> & gs_jets = ghosted_seq.
jets();
574 const double tolerance = 1e-11;
579 valarray<double> our_areas(_history.size());
582 valarray<PseudoJet> our_area_4vectors(_history.size());
583 our_area_4vectors =
PseudoJet(0.0,0.0,0.0,0.0);
585 for (
unsigned i = 0; i < gs_history.size(); i++) {
587 unsigned gs_hist_index = gs_unique_hist_order[i];
588 if (gs_hist_index < ghosted_seq.
n_particles())
continue;
590 int parent1 = gs_hist.parent1;
593 if (parent2 == BeamJet) {
596 gs_jets[gs_history[parent1].jetp_index];
597 double area_local = ghosted_seq.
area(jet);
602 _ghost_jets.push_back(GhostJet(jet,area_local));
603 if (abs(jet.
rap()) < _safe_rap_for_area) {
604 _non_jet_area += area_local;
605 _non_jet_area2 += area_local*area_local;
606 _non_jet_number += 1;
613 while (++j <
int(_history.size())) {
614 hist_index = unique_hist_order[j];
615 if (hist_index >= _initial_n)
break;}
618 if (j >=
int(_history.size()))
throw Error(
"ClusterSequenceActiveArea: overran reference array in diB matching");
622 int refjet_index = _history[_history[hist_index].parent1].jetp_index;
623 assert(refjet_index >= 0 && refjet_index <
int(_jets.size()));
624 const PseudoJet & refjet = _jets[refjet_index];
634 _throw_unless_jets_have_same_perp_or_E(jet, refjet, tolerance,
638 our_areas[hist_index] = area_local;
639 our_area_4vectors[hist_index] = ext_area;
645 our_areas[_history[hist_index].parent1] = area_local;
646 our_area_4vectors[_history[hist_index].parent1] = ext_area;
654 while (++j <
int(_history.size())) {
655 hist_index = unique_hist_order[j];
656 if (hist_index >= _initial_n)
break;}
659 if (j >=
int(_history.size()))
throw Error(
"ClusterSequenceActiveArea: overran reference array in dij matching");
664 if (_history[hist_index].parent2 == BeamJet)
throw Error(
"ClusterSequenceActiveArea: could not match clustering sequences (encountered dij matched with diB)");
671 const PseudoJet & refjet = _jets[_history[hist_index].jetp_index];
674 _throw_unless_jets_have_same_perp_or_E(jet, refjet, tolerance,
679 double area_local = ghosted_seq.
area(jet);
680 our_areas[hist_index] += area_local;
699 _jet_def.recombiner()->plus_equal(our_area_4vectors[hist_index], ext_area);
707 const PseudoJet & jet1 = gs_jets[gs_history[parent1].jetp_index];
708 int our_parent1 = _history[hist_index].parent1;
709 our_areas[our_parent1] = ghosted_seq.
area(jet1);
710 our_area_4vectors[our_parent1] = ghosted_seq.
area_4vector(jet1);
712 const PseudoJet & jet2 = gs_jets[gs_history[parent2].jetp_index];
713 int our_parent2 = _history[hist_index].parent2;
714 our_areas[our_parent2] = ghosted_seq.
area(jet2);
715 our_area_4vectors[our_parent2] = ghosted_seq.
area_4vector(jet2);
723 for (
unsigned iu = 0; iu < unclust.size(); iu++) {
725 double area_local = ghosted_seq.
area(unclust[iu]);
726 _unclustered_ghosts.push_back(GhostJet(unclust[iu],area_local));
738 for (
unsigned int area_index = 0; area_index<our_areas.size(); area_index++){
739 _average_area[area_index] += our_areas[area_index];
740 _average_area2[area_index] += our_areas[area_index]*our_areas[area_index];
746 for (
unsigned i = 0; i < our_area_4vectors.size(); i++) {
747 _jet_def.recombiner()->plus_equal(_average_area_4vector[i],
748 our_area_4vectors[i]);
756 void ClusterSequenceActiveArea::_throw_unless_jets_have_same_perp_or_E(
765 && abs(jet.E()-refjet.E()) > tolerance*max(jet.E(),refjet.E())) {
767 ostr <<
"Could not match clustering sequence for an inclusive/exclusive jet when reconstructing areas. See FAQ for possible explanations." << endl;
769 << refjet.px() <<
" " 770 << refjet.py() <<
" " 771 << refjet.pz() <<
" " 772 << refjet.E() << endl;
779 ostr <<
" NB: some particles have pt too low wrt ghosts -- this may be the cause" << endl;}
782 throw Error(ostr.str());
786 FASTJET_END_NAMESPACE
virtual bool is_pure_ghost(const PseudoJet &jet) const override
true if a jet is made exclusively of ghosts
std::vector< int > unique_history_order() const
routine that returns an order in which to read the history such that clusterings that lead to identic...
mean_pt_strategies
enum providing a variety of tentative strategies for estimating the background (e.g.
const std::vector< history_element > & history() const
allow the user to access the raw internal history.
int jetp_index
index in _history where the current jet is recombined with another jet to form its child...
virtual double area(const PseudoJet &jet) const override
returns the area of a jet
int parent2
index in _history where first parent of this jet was created (InexistentParent if this jet is an orig...
bool applies_jet_by_jet() const
returns true if this can be applied jet by jet
double dij
index in the _jets vector where we will find the
bool pass(const PseudoJet &jet) const
return true if the jet passes the selection
Strategy strategy_used() const
return the enum value of the strategy used to cluster the event
base class corresponding to errors that can be thrown by FastJet
Like ClusterSequence with computation of the active jet area with the addition of explicit ghosts...
double perp2() const
returns the squared transverse momentum
bool has_dangerous_particles() const
returns true if there are any particles whose transverse momentum if so low that there's a risk of th...
double rap() const
returns the rapidity or some large value when the rapidity is infinite
Class that encodes information about cuts and other selection criteria that can be applied to PseudoJ...
a single element in the clustering history
virtual PseudoJet area_4vector(const PseudoJet &jet) const override
returns a four vector corresponding to the sum (E-scheme) of the ghost four-vectors composing the jet...
Parameters to configure the computation of jet areas using ghosts.
Class to contain pseudojets, including minimal information of use to jet-clustering routines...
const std::vector< PseudoJet > & jets() const
allow the user to access the internally stored _jets() array, which contains both the initial particl...
class that is intended to hold a full definition of the jet clusterer
unsigned int n_particles() const
returns the number of particles that were provided to the clustering algorithm (helps the user find t...
std::vector< PseudoJet > unclustered_particles() const
return the set of particles that have not been clustered.