1 | //STARTHEADER
|
---|
2 | // $Id: MassDropTagger.cc 1332 2013-11-20 20:52:59Z pavel $
|
---|
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 | #include <fastjet/tools/MassDropTagger.hh>
|
---|
30 | #include <fastjet/ClusterSequence.hh>
|
---|
31 | #include <sstream>
|
---|
32 |
|
---|
33 | FASTJET_BEGIN_NAMESPACE
|
---|
34 |
|
---|
35 | LimitedWarning MassDropTagger::_warnings_nonca;
|
---|
36 |
|
---|
37 | using namespace std;
|
---|
38 |
|
---|
39 | //----------------------------------------------------------------------
|
---|
40 | // MassDropTagger class implementation
|
---|
41 | //----------------------------------------------------------------------
|
---|
42 |
|
---|
43 | //------------------------------------------------------------------------
|
---|
44 | // description of the tagger
|
---|
45 | string MassDropTagger::description() const{
|
---|
46 | ostringstream oss;
|
---|
47 | oss << "MassDropTagger with mu=" << _mu << " and ycut=" << _ycut;
|
---|
48 | return oss.str();
|
---|
49 | }
|
---|
50 |
|
---|
51 | //------------------------------------------------------------------------
|
---|
52 | // returns the tagged PseudoJet if successful, 0 otherwise
|
---|
53 | // - jet the PseudoJet to tag
|
---|
54 | PseudoJet MassDropTagger::result(const PseudoJet & jet) const{
|
---|
55 | PseudoJet j = jet;
|
---|
56 |
|
---|
57 | // issue a warning if the jet is not obtained through a C/A
|
---|
58 | // clustering
|
---|
59 | if ((! j.has_associated_cluster_sequence()) ||
|
---|
60 | (j.validated_cs()->jet_def().jet_algorithm() != cambridge_algorithm))
|
---|
61 | _warnings_nonca.warn("MassDropTagger should only be applied on jets from a Cambridge/Aachen clustering; use it with other algorithms at your own risk.");
|
---|
62 |
|
---|
63 |
|
---|
64 | PseudoJet j1, j2;
|
---|
65 | bool had_parents;
|
---|
66 |
|
---|
67 | // we just ask that we can "walk" in the cluster sequence.
|
---|
68 | // appropriate errors will be thrown automatically if this is not
|
---|
69 | // the case
|
---|
70 | while ((had_parents = j.has_parents(j1,j2))) {
|
---|
71 | // make parent1 the more massive jet
|
---|
72 | if (j1.m2() < j2.m2()) std::swap(j1,j2);
|
---|
73 |
|
---|
74 | // if we pass the conditions on the mass drop and its degree of
|
---|
75 | // asymmetry (kt_dist/m^2 > rtycut [where kt_dist/m^2 \sim
|
---|
76 | // z/(1-z)), then we've found something interesting, so exit the
|
---|
77 | // loop
|
---|
78 | if ( (j1.m2() < _mu*_mu*j.m2()) && (j1.kt_distance(j2) > _ycut*j.m2()) )
|
---|
79 | break;
|
---|
80 | else
|
---|
81 | j = j1;
|
---|
82 | }
|
---|
83 |
|
---|
84 | if (!had_parents)
|
---|
85 | // no Higgs found, return an empty PseudoJet
|
---|
86 | return PseudoJet();
|
---|
87 |
|
---|
88 | // create the result and its structure
|
---|
89 | PseudoJet result_local = j;
|
---|
90 | MassDropTaggerStructure * s = new MassDropTaggerStructure(result_local);
|
---|
91 | // s->_original_jet = jet;
|
---|
92 | s->_mu = (j.m2()!=0.0) ? sqrt(j1.m2()/j.m2()) : 0.0;
|
---|
93 | s->_y = (j.m2()!=0.0) ? j1.kt_distance(j2)/j.m2() : 0.0;
|
---|
94 |
|
---|
95 | result_local.set_structure_shared_ptr(SharedPtr<PseudoJetStructureBase>(s));
|
---|
96 |
|
---|
97 | return result_local;
|
---|
98 | }
|
---|
99 |
|
---|
100 | FASTJET_END_NAMESPACE
|
---|
101 |
|
---|