CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Protected Member Functions | Private Attributes
contrib::CMSBoostedTauSeedingAlgorithm Class Reference

#include <CMSBoostedTauSeedingAlgorithm.h>

Inheritance diagram for contrib::CMSBoostedTauSeedingAlgorithm:

Public Types

typedef CMSBoostedTauSeedingAlgorithmStructure StructureType
 

Public Member Functions

 CMSBoostedTauSeedingAlgorithm (double ptMin, double muMin, double muMax, double yMin, double yMax, double dRMin, double dRMax, int maxDepth, int verbosity=0)
 
std::string description () const override
 
PseudoJet result (const PseudoJet &jet) const override
 
 ~CMSBoostedTauSeedingAlgorithm () override
 

Protected Member Functions

void dumpSubJetStructure (const fastjet::PseudoJet &jet, int depth, int maxDepth, const std::string &depth_and_idx_string) const
 
std::pair< PseudoJet, PseudoJet > findSubjets (const PseudoJet &jet, int depth, bool &subjetsFound) const
 

Private Attributes

double dRMax_
 the max value of the dR parameter More...
 
double dRMin_
 the min value of the dR parameter More...
 
int maxDepth_
 the max depth for descending into clustering sequence More...
 
double muMax_
 the max value of the mass-drop parameter More...
 
double muMin_
 the min value of the mass-drop parameter More...
 
double ptMin_
 minimum sub-jet pt More...
 
int verbosity_
 flag to enable/disable debug output More...
 
double yMax_
 the max value of the asymmetry parameter More...
 
double yMin_
 the min value of the asymmetry parameter More...
 

Detailed Description

This class implements the CMS boosted tau algorithm

Definition at line 52 of file CMSBoostedTauSeedingAlgorithm.h.

Member Typedef Documentation

◆ StructureType

Definition at line 74 of file CMSBoostedTauSeedingAlgorithm.h.

Constructor & Destructor Documentation

◆ CMSBoostedTauSeedingAlgorithm()

contrib::CMSBoostedTauSeedingAlgorithm::CMSBoostedTauSeedingAlgorithm ( double  ptMin,
double  muMin,
double  muMax,
double  yMin,
double  yMax,
double  dRMin,
double  dRMax,
int  maxDepth,
int  verbosity = 0 
)

Definition at line 30 of file CMSBoostedTauSeedingAlgorithm.cc.

39  : ptMin_(iminPt),
40  muMin_(iminMassDrop),
41  muMax_(imaxMassDrop),
42  yMin_(iminY),
43  yMax_(imaxY),
44  dRMin_(iminDeltaR),
45  dRMax_(imaxDeltaR),

◆ ~CMSBoostedTauSeedingAlgorithm()

contrib::CMSBoostedTauSeedingAlgorithm::~CMSBoostedTauSeedingAlgorithm ( )
inlineoverride

Definition at line 66 of file CMSBoostedTauSeedingAlgorithm.h.

66 {}

Member Function Documentation

◆ description()

std::string contrib::CMSBoostedTauSeedingAlgorithm::description ( ) const
override

Definition at line 51 of file CMSBoostedTauSeedingAlgorithm.cc.

51  {
52  std::ostringstream oss;
53  oss << "CMSBoostedTauSeedingAlgorithm algorithm";
54  return oss.str();
55  }

◆ dumpSubJetStructure()

void contrib::CMSBoostedTauSeedingAlgorithm::dumpSubJetStructure ( const fastjet::PseudoJet &  jet,
int  depth,
int  maxDepth,
const std::string &  depth_and_idx_string 
) const
protected

Definition at line 58 of file CMSBoostedTauSeedingAlgorithm.cc.

59  {
60  if (maxDepth != -1 && depth > maxDepth)
61  return;
62  fastjet::PseudoJet subjet1, subjet2;
63  bool hasSubjets = jet.has_parents(subjet1, subjet2);
64  if (!hasSubjets)
65  return;
66  std::string depth_and_idx_string_subjet1 = depth_and_idx_string;
67  if (depth_and_idx_string_subjet1.length() > 0)
68  depth_and_idx_string_subjet1.append(".");
69  depth_and_idx_string_subjet1.append("0");
70  for (int iSpace = 0; iSpace < depth; ++iSpace) {
71  std::cout << " ";
72  }
73  std::cout << " jetConstituent #" << depth_and_idx_string_subjet1 << " (depth = " << depth
74  << "): Pt = " << subjet1.pt() << ","
75  << " eta = " << subjet1.eta() << ", phi = " << subjet1.phi() << ", mass = " << subjet1.m()
76  << " (constituents = " << subjet1.constituents().size() << ")" << std::endl;
77  dumpSubJetStructure(subjet1, depth + 1, maxDepth, depth_and_idx_string_subjet1);
78  std::string depth_and_idx_string_subjet2 = depth_and_idx_string;
79  if (depth_and_idx_string_subjet2.length() > 0)
80  depth_and_idx_string_subjet2.append(".");
81  depth_and_idx_string_subjet2.append("1");
82  for (int iSpace = 0; iSpace < depth; ++iSpace) {
83  std::cout << " ";
84  }
85  std::cout << " jetConstituent #" << depth_and_idx_string_subjet2 << " (depth = " << depth
86  << "): Pt = " << subjet2.pt() << ","
87  << " eta = " << subjet2.eta() << ", phi = " << subjet2.phi() << ", mass = " << subjet2.m()
88  << " (constituents = " << subjet2.constituents().size() << ")" << std::endl;
89  dumpSubJetStructure(subjet2, depth + 1, maxDepth, depth_and_idx_string_subjet2);
90  for (int iSpace = 0; iSpace < depth; ++iSpace) {
91  std::cout << " ";
92  }
93  double dR = subjet1.delta_R(subjet2);
94  std::cout << " (mass-drop @ " << depth_and_idx_string << " = " << std::max(subjet1.m(), subjet2.m()) / jet.m()
95  << ", dR = " << dR << ")" << std::endl;
96  }

References gather_cfg::cout, LEDCalibrationChannels::depth, HGC3DClusterGenMatchSelector_cfi::dR, dumpSubJetStructure(), metsig::jet, SiStripPI::max, HLT_FULL_cff::maxDepth, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by dumpSubJetStructure(), and result().

◆ findSubjets()

std::pair< PseudoJet, PseudoJet > contrib::CMSBoostedTauSeedingAlgorithm::findSubjets ( const PseudoJet &  jet,
int  depth,
bool &  subjetsFound 
) const
protected

Definition at line 98 of file CMSBoostedTauSeedingAlgorithm.cc.

99  {
100  const float etaMax_ = 3.0; // cut on eta subjets <3. Anyway taus are cut at 2.3
101  if (verbosity_ >= 2) {
102  std::cout << "<CMSBoostedTauSeedingAlgorithm::findSubjets>:" << std::endl;
103  std::cout << " jet: Pt = " << jet.pt() << ", eta = " << jet.eta() << ", phi = " << jet.phi()
104  << ", mass = " << jet.m() << std::endl;
105  }
106 
107  PseudoJet subjet1, subjet2;
108  bool hasSubjets = jet.has_parents(subjet1, subjet2);
109  if (hasSubjets && (maxDepth_ == -1 || depth <= maxDepth_)) {
110  // make subjet1 the more massive jet
111  if (subjet1.m2() < subjet2.m2()) {
112  std::swap(subjet1, subjet2);
113  }
114  double dR = subjet1.delta_R(subjet2);
115  double kT = subjet1.kt_distance(subjet2);
116  double mu = (jet.m() > 0.) ? sqrt(std::max(subjet1.m2(), subjet2.m2()) / jet.m2()) : -1.;
117  // check if subjets pass selection required for seeding boosted tau reconstruction
118  if (subjet1.pt() > ptMin_ && fabs(subjet1.eta()) < etaMax_ && subjet2.pt() > ptMin_ &&
119  fabs(subjet2.eta()) < etaMax_ && dR > dRMin_ && dR < dRMax_ && mu > muMin_ && mu < muMax_ &&
120  kT < (yMax_ * jet.m2()) && kT > (yMin_ * jet.m2())) {
121  subjetsFound = true;
122  return std::make_pair(subjet1, subjet2);
123  } else if (subjet1.pt() > ptMin_) {
124  return findSubjets(subjet1, depth + 1, subjetsFound);
125  } else if (subjet2.pt() > ptMin_) {
126  return findSubjets(subjet2, depth + 1, subjetsFound);
127  }
128  }
129  subjetsFound = false;
130  PseudoJet dummy_subjet1, dummy_subjet2;
131  return std::make_pair(dummy_subjet1, dummy_subjet2);
132  }

References gather_cfg::cout, LEDCalibrationChannels::depth, HGC3DClusterGenMatchSelector_cfi::dR, dRMin_, metsig::jet, SiStripPI::max, maxDepth_, amptDefaultParameters_cff::mu, muMax_, muMin_, ptMin_, mathSSE::sqrt(), std::swap(), verbosity_, yMax_, and yMin_.

Referenced by result().

◆ result()

PseudoJet contrib::CMSBoostedTauSeedingAlgorithm::result ( const PseudoJet &  jet) const
override

Definition at line 134 of file CMSBoostedTauSeedingAlgorithm.cc.

134  {
135  if (verbosity_ >= 1) {
136  std::cout << "<CMSBoostedTauSeedingAlgorithm::findSubjets>:" << std::endl;
137  std::cout << " jet: Pt = " << jet.pt() << ", eta = " << jet.eta() << ", phi = " << jet.phi()
138  << ", mass = " << jet.m() << std::endl;
139  }
140 
141  if (verbosity_ >= 2) {
143  }
144 
145  bool subjetsFound = false;
146  std::pair<PseudoJet, PseudoJet> subjets = findSubjets(jet, 0, subjetsFound);
147  if (subjetsFound) {
148  // fill structure for returning result
149  PseudoJet subjet1 = subjets.first;
150  PseudoJet subjet2 = subjets.second;
151  if (verbosity_ >= 1) {
152  std::cout << "before recombination:" << std::endl;
153  std::cout << " subjet #1: Pt = " << subjet1.pt() << ", eta = " << subjet1.eta() << ", phi = " << subjet1.phi()
154  << ", mass = " << subjet1.m() << std::endl;
155  std::cout << " subjet #2: Pt = " << subjet2.pt() << ", eta = " << subjet2.eta() << ", phi = " << subjet2.phi()
156  << ", mass = " << subjet2.m() << std::endl;
157  }
158 
159  const JetDefinition::Recombiner* rec = jet.associated_cluster_sequence()->jet_def().recombiner();
160  PseudoJet result_local = join(subjet1, subjet2, *rec);
161  if (verbosity_ >= 1) {
162  std::cout << "after recombination:" << std::endl;
163  std::vector<fastjet::PseudoJet> subjets = result_local.pieces();
164  int idx_subjet = 0;
165  for (std::vector<fastjet::PseudoJet>::const_iterator subjet = subjets.begin(); subjet != subjets.end();
166  ++subjet) {
167  std::cout << " subjet #" << idx_subjet << ": Pt = " << subjet->pt() << ", eta = " << subjet->eta()
168  << ", phi = " << subjet->phi() << ", mass = " << subjet->m()
169  << " (#constituents = " << subjet->constituents().size() << ")" << std::endl;
170  std::vector<fastjet::PseudoJet> constituents = subjet->constituents();
171  int idx_constituent = 0;
172  for (std::vector<fastjet::PseudoJet>::const_iterator constituent = constituents.begin();
173  constituent != constituents.end();
174  ++constituent) {
175  if (constituent->pt() < 1.e-3)
176  continue; // CV: skip ghosts
177  std::cout << " constituent #" << idx_constituent << ": Pt = " << constituent->pt()
178  << ", eta = " << constituent->eta() << ", phi = " << constituent->phi() << ","
179  << " mass = " << constituent->m() << std::endl;
180  ++idx_constituent;
181  }
182  ++idx_subjet;
183  }
184  }
185 
186  CMSBoostedTauSeedingAlgorithmStructure* s = new CMSBoostedTauSeedingAlgorithmStructure(result_local);
187  //s->_original_jet = jet;
188  s->_mu = (jet.m2() > 0.) ? sqrt(std::max(subjet1.m2(), subjet2.m2()) / jet.m2()) : 0.;
189  s->_y = (jet.m2() > 0.) ? subjet1.kt_distance(subjet2) / jet.m2() : 0.;
190  s->_dR = subjet1.delta_R(subjet2);
191  s->_pt = subjet2.pt();
192 
193  result_local.set_structure_shared_ptr(SharedPtr<PseudoJetStructureBase>(s));
194 
195  return result_local;
196  } else {
197  // no subjets for seeding boosted tau reconstruction found, return an empty PseudoJet
198  if (verbosity_ >= 1) {
199  std::cout << "No subjets found." << std::endl;
200  }
201  return PseudoJet();
202  }
203  }

References gather_cfg::cout, dumpSubJetStructure(), findSubjets(), metsig::jet, join(), SiStripPI::max, maxDepth_, alignCSCRings::s, mathSSE::sqrt(), and verbosity_.

Member Data Documentation

◆ dRMax_

double contrib::CMSBoostedTauSeedingAlgorithm::dRMax_
private

the max value of the dR parameter

Definition at line 90 of file CMSBoostedTauSeedingAlgorithm.h.

◆ dRMin_

double contrib::CMSBoostedTauSeedingAlgorithm::dRMin_
private

the min value of the dR parameter

Definition at line 89 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets().

◆ maxDepth_

int contrib::CMSBoostedTauSeedingAlgorithm::maxDepth_
private

the max depth for descending into clustering sequence

Definition at line 91 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets(), and result().

◆ muMax_

double contrib::CMSBoostedTauSeedingAlgorithm::muMax_
private

the max value of the mass-drop parameter

Definition at line 86 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets().

◆ muMin_

double contrib::CMSBoostedTauSeedingAlgorithm::muMin_
private

the min value of the mass-drop parameter

Definition at line 85 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets().

◆ ptMin_

double contrib::CMSBoostedTauSeedingAlgorithm::ptMin_
private

minimum sub-jet pt

Definition at line 84 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets().

◆ verbosity_

int contrib::CMSBoostedTauSeedingAlgorithm::verbosity_
private

flag to enable/disable debug output

Definition at line 93 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets(), and result().

◆ yMax_

double contrib::CMSBoostedTauSeedingAlgorithm::yMax_
private

the max value of the asymmetry parameter

Definition at line 88 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets().

◆ yMin_

double contrib::CMSBoostedTauSeedingAlgorithm::yMin_
private

the min value of the asymmetry parameter

Definition at line 87 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets().

HIPAlignmentAlgorithm_cfi.verbosity
verbosity
Definition: HIPAlignmentAlgorithm_cfi.py:7
contrib::CMSBoostedTauSeedingAlgorithm::muMin_
double muMin_
the min value of the mass-drop parameter
Definition: CMSBoostedTauSeedingAlgorithm.h:85
contrib::CMSBoostedTauSeedingAlgorithm::muMax_
double muMax_
the max value of the mass-drop parameter
Definition: CMSBoostedTauSeedingAlgorithm.h:86
amptDefaultParameters_cff.mu
mu
Definition: amptDefaultParameters_cff.py:16
gather_cfg.cout
cout
Definition: gather_cfg.py:144
join
static std::string join(char **cmd)
Definition: RemoteFile.cc:17
contrib::CMSBoostedTauSeedingAlgorithm::ptMin_
double ptMin_
minimum sub-jet pt
Definition: CMSBoostedTauSeedingAlgorithm.h:84
contrib::CMSBoostedTauSeedingAlgorithm::yMax_
double yMax_
the max value of the asymmetry parameter
Definition: CMSBoostedTauSeedingAlgorithm.h:88
contrib::CMSBoostedTauSeedingAlgorithm::verbosity_
int verbosity_
flag to enable/disable debug output
Definition: CMSBoostedTauSeedingAlgorithm.h:93
HLT_FULL_cff.maxDepth
maxDepth
Definition: HLT_FULL_cff.py:8669
alignCSCRings.s
s
Definition: alignCSCRings.py:92
std::swap
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
Definition: DataFrameContainer.h:209
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
contrib::CMSBoostedTauSeedingAlgorithm::dRMin_
double dRMin_
the min value of the dR parameter
Definition: CMSBoostedTauSeedingAlgorithm.h:89
LEDCalibrationChannels.depth
depth
Definition: LEDCalibrationChannels.py:65
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
contrib::CMSBoostedTauSeedingAlgorithm::dumpSubJetStructure
void dumpSubJetStructure(const fastjet::PseudoJet &jet, int depth, int maxDepth, const std::string &depth_and_idx_string) const
Definition: CMSBoostedTauSeedingAlgorithm.cc:58
contrib::CMSBoostedTauSeedingAlgorithm::findSubjets
std::pair< PseudoJet, PseudoJet > findSubjets(const PseudoJet &jet, int depth, bool &subjetsFound) const
Definition: CMSBoostedTauSeedingAlgorithm.cc:98
contrib::CMSBoostedTauSeedingAlgorithm::dRMax_
double dRMax_
the max value of the dR parameter
Definition: CMSBoostedTauSeedingAlgorithm.h:90
metsig::jet
Definition: SignAlgoResolutions.h:47
HGC3DClusterGenMatchSelector_cfi.dR
dR
Definition: HGC3DClusterGenMatchSelector_cfi.py:7
contrib::CMSBoostedTauSeedingAlgorithm::maxDepth_
int maxDepth_
the max depth for descending into clustering sequence
Definition: CMSBoostedTauSeedingAlgorithm.h:91
contrib::CMSBoostedTauSeedingAlgorithm::yMin_
double yMin_
the min value of the asymmetry parameter
Definition: CMSBoostedTauSeedingAlgorithm.h:87