CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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)
 
virtual std::string description () const
 
virtual PseudoJet result (const PseudoJet &jet) const
 
virtual ~CMSBoostedTauSeedingAlgorithm ()
 

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 53 of file CMSBoostedTauSeedingAlgorithm.h.

Member Typedef Documentation

Definition at line 68 of file CMSBoostedTauSeedingAlgorithm.h.

Constructor & Destructor Documentation

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 29 of file CMSBoostedTauSeedingAlgorithm.cc.

35  : ptMin_(iminPt),
36  muMin_(iminMassDrop), muMax_(imaxMassDrop),
37  yMin_(iminY), yMax_(imaxY),
38  dRMin_(iminDeltaR), dRMax_(imaxDeltaR),
39  maxDepth_(maxDepth),
41  {}
double dRMax_
the max value of the dR parameter
int maxDepth_
the max depth for descending into clustering sequence
double yMax_
the max value of the asymmetry parameter
double yMin_
the min value of the asymmetry parameter
double muMax_
the max value of the mass-drop parameter
int verbosity_
flag to enable/disable debug output
double dRMin_
the min value of the dR parameter
double muMin_
the min value of the mass-drop parameter
tuple verbosity
Definition: mvaPFMET_cff.py:80
virtual contrib::CMSBoostedTauSeedingAlgorithm::~CMSBoostedTauSeedingAlgorithm ( )
inlinevirtual

Definition at line 60 of file CMSBoostedTauSeedingAlgorithm.h.

60 {}

Member Function Documentation

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

Definition at line 45 of file CMSBoostedTauSeedingAlgorithm.cc.

45  {
46  std::ostringstream oss;
47  oss << "CMSBoostedTauSeedingAlgorithm algorithm";
48  return oss.str();
49  }
void contrib::CMSBoostedTauSeedingAlgorithm::dumpSubJetStructure ( const fastjet::PseudoJet &  jet,
int  depth,
int  maxDepth,
const std::string &  depth_and_idx_string 
) const
protected

Definition at line 53 of file CMSBoostedTauSeedingAlgorithm.cc.

References gather_cfg::cout, PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, bookConverter::max, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by result().

54  {
55  if ( maxDepth != -1 && depth > maxDepth ) return;
56  fastjet::PseudoJet subjet1, subjet2;
57  bool hasSubjets = jet.has_parents(subjet1, subjet2);
58  if ( !hasSubjets ) return;
59  std::string depth_and_idx_string_subjet1 = depth_and_idx_string;
60  if ( depth_and_idx_string_subjet1.length() > 0 ) depth_and_idx_string_subjet1.append(".");
61  depth_and_idx_string_subjet1.append("0");
62  for ( int iSpace = 0; iSpace < depth; ++iSpace ) {
63  std::cout << " ";
64  }
65  std::cout << " jetConstituent #" << depth_and_idx_string_subjet1 << " (depth = " << depth << "): Pt = " << subjet1.pt() << ","
66  << " eta = " << subjet1.eta() << ", phi = " << subjet1.phi() << ", mass = " << subjet1.m()
67  << " (constituents = " << subjet1.constituents().size() << ")" << std::endl;
68  dumpSubJetStructure(subjet1, depth + 1, maxDepth, depth_and_idx_string_subjet1);
69  std::string depth_and_idx_string_subjet2 = depth_and_idx_string;
70  if ( depth_and_idx_string_subjet2.length() > 0 ) depth_and_idx_string_subjet2.append(".");
71  depth_and_idx_string_subjet2.append("1");
72  for ( int iSpace = 0; iSpace < depth; ++iSpace ) {
73  std::cout << " ";
74  }
75  std::cout << " jetConstituent #" << depth_and_idx_string_subjet2 << " (depth = " << depth << "): Pt = " << subjet2.pt() << ","
76  << " eta = " << subjet2.eta() << ", phi = " << subjet2.phi() << ", mass = " << subjet2.m()
77  << " (constituents = " << subjet2.constituents().size() << ")" << std::endl;
78  dumpSubJetStructure(subjet2, depth + 1, maxDepth, depth_and_idx_string_subjet2);
79  for ( int iSpace = 0; iSpace < depth; ++iSpace ) {
80  std::cout << " ";
81  }
82  double dR = subjet1.delta_R(subjet2);
83  std::cout << " (mass-drop @ " << depth_and_idx_string << " = " << std::max(subjet1.m(), subjet2.m())/jet.m() << ", dR = " << dR << ")" << std::endl;
84  }
void dumpSubJetStructure(const fastjet::PseudoJet &jet, int depth, int maxDepth, const std::string &depth_and_idx_string) const
tuple cout
Definition: gather_cfg.py:121
std::pair< PseudoJet, PseudoJet > contrib::CMSBoostedTauSeedingAlgorithm::findSubjets ( const PseudoJet &  jet,
int  depth,
bool &  subjetsFound 
) const
protected

Definition at line 86 of file CMSBoostedTauSeedingAlgorithm.cc.

References gather_cfg::cout, PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, dRMin_, bookConverter::max, maxDepth_, RPCpg::mu, muMax_, muMin_, ptMin_, mathSSE::sqrt(), std::swap(), verbosity_, yMax_, and yMin_.

Referenced by result().

87  {
88  if ( verbosity_ >= 2 ) {
89  std::cout << "<CMSBoostedTauSeedingAlgorithm::findSubjets>:" << std::endl;
90  std::cout << " jet: Pt = " << jet.pt() << ", eta = " << jet.eta() << ", phi = " << jet.phi() << ", mass = " << jet.m() << std::endl;
91  }
92 
93  PseudoJet subjet1, subjet2;
94  bool hasSubjets = jet.has_parents(subjet1, subjet2);
95  if ( hasSubjets && (maxDepth_ == -1 || depth <= maxDepth_) ) {
96  // make subjet1 the more massive jet
97  if ( subjet1.m2() < subjet2.m2() ) {
98  std::swap(subjet1, subjet2);
99  }
100  double dR = subjet1.delta_R(subjet2);
101  double kT = subjet1.kt_distance(subjet2);
102  double mu = ( jet.m() > 0. ) ?
103  sqrt(std::max(subjet1.m2(), subjet2.m2())/jet.m2()) : -1.;
104  // check if subjets pass selection required for seeding boosted tau reconstruction
105  if ( subjet1.pt() > ptMin_ && subjet2.pt() > ptMin_ && dR > dRMin_ && dR < dRMax_ && mu > muMin_ && mu < muMax_ && kT < (yMax_*jet.m2()) && kT > (yMin_*jet.m2()) ) {
106  subjetsFound = true;
107  return std::make_pair(subjet1, subjet2);
108  } else if ( subjet1.pt() > ptMin_ ) {
109  return findSubjets(subjet1, depth + 1, subjetsFound);
110  } else if ( subjet2.pt() > ptMin_ ) {
111  return findSubjets(subjet2, depth + 1, subjetsFound);
112  }
113  }
114  subjetsFound = false;
115  PseudoJet dummy_subjet1, dummy_subjet2;
116  return std::make_pair(dummy_subjet1, dummy_subjet2);
117  }
int maxDepth_
the max depth for descending into clustering sequence
double yMax_
the max value of the asymmetry parameter
void swap(edm::DataFrameContainer &lhs, edm::DataFrameContainer &rhs)
T sqrt(T t)
Definition: SSEVec.h:48
double yMin_
the min value of the asymmetry parameter
double muMax_
the max value of the mass-drop parameter
const int mu
Definition: Constants.h:22
std::pair< PseudoJet, PseudoJet > findSubjets(const PseudoJet &jet, int depth, bool &subjetsFound) const
int verbosity_
flag to enable/disable debug output
double dRMin_
the min value of the dR parameter
tuple cout
Definition: gather_cfg.py:121
double muMin_
the min value of the mass-drop parameter
PseudoJet contrib::CMSBoostedTauSeedingAlgorithm::result ( const PseudoJet &  jet) const
virtual

Definition at line 119 of file CMSBoostedTauSeedingAlgorithm.cc.

References contrib::CMSBoostedTauSeedingAlgorithmStructure::_dR, contrib::CMSBoostedTauSeedingAlgorithmStructure::_mu, contrib::CMSBoostedTauSeedingAlgorithmStructure::_pt, contrib::CMSBoostedTauSeedingAlgorithmStructure::_y, gather_cfg::cout, dumpSubJetStructure(), findSubjets(), join(), bookConverter::max, maxDepth_, alignCSCRings::s, mathSSE::sqrt(), and verbosity_.

120  {
121  if ( verbosity_ >= 1 ) {
122  std::cout << "<CMSBoostedTauSeedingAlgorithm::findSubjets>:" << std::endl;
123  std::cout << " jet: Pt = " << jet.pt() << ", eta = " << jet.eta() << ", phi = " << jet.phi() << ", mass = " << jet.m() << std::endl;
124  }
125 
126  if ( verbosity_ >= 2 ) {
128  }
129 
130  bool subjetsFound = false;
131  std::pair<PseudoJet, PseudoJet> subjets = findSubjets(jet, 0, subjetsFound);
132  if ( subjetsFound ) {
133  // fill structure for returning result
134  PseudoJet subjet1 = subjets.first;
135  PseudoJet subjet2 = subjets.second;
136  if ( verbosity_ >= 1 ) {
137  std::cout << "before recombination:" << std::endl;
138  std::cout << " subjet #1: Pt = " << subjet1.pt() << ", eta = " << subjet1.eta() << ", phi = " << subjet1.phi() << ", mass = " << subjet1.m() << std::endl;
139  std::cout << " subjet #2: Pt = " << subjet2.pt() << ", eta = " << subjet2.eta() << ", phi = " << subjet2.phi() << ", mass = " << subjet2.m() << std::endl;
140  }
141 
142  const JetDefinition::Recombiner* rec = jet.associated_cluster_sequence()->jet_def().recombiner();
143  PseudoJet result_local = join(subjet1, subjet2, *rec);
144  if ( verbosity_ >= 1 ) {
145  std::cout << "after recombination:" << std::endl;
146  std::vector<fastjet::PseudoJet> subjets = result_local.pieces();
147  int idx_subjet = 0;
148  for ( std::vector<fastjet::PseudoJet>::const_iterator subjet = subjets.begin();
149  subjet != subjets.end(); ++subjet ) {
150  std::cout << " subjet #" << idx_subjet << ": Pt = " << subjet->pt() << ", eta = " << subjet->eta() << ", phi = " << subjet->phi() << ", mass = " << subjet->m()
151  << " (#constituents = " << subjet->constituents().size() << ")" << std::endl;
152  std::vector<fastjet::PseudoJet> constituents = subjet->constituents();
153  int idx_constituent = 0;
154  for ( std::vector<fastjet::PseudoJet>::const_iterator constituent = constituents.begin();
155  constituent != constituents.end(); ++constituent ) {
156  if ( constituent->pt() < 1.e-3 ) continue; // CV: skip ghosts
157  std::cout << " constituent #" << idx_constituent << ": Pt = " << constituent->pt() << ", eta = " << constituent->eta() << ", phi = " << constituent->phi() << ","
158  << " mass = " << constituent->m() << std::endl;
159  ++idx_constituent;
160  }
161  ++idx_subjet;
162  }
163  }
164 
165  CMSBoostedTauSeedingAlgorithmStructure* s = new CMSBoostedTauSeedingAlgorithmStructure(result_local);
166  //s->_original_jet = jet;
167  s->_mu = ( jet.m2() > 0. ) ? sqrt(std::max(subjet1.m2(), subjet2.m2())/jet.m2()) : 0.;
168  s->_y = ( jet.m2() > 0. ) ? subjet1.kt_distance(subjet2)/jet.m2() : 0.;
169  s->_dR = subjet1.delta_R(subjet2);
170  s->_pt = subjet2.pt();
171 
172  result_local.set_structure_shared_ptr(SharedPtr<PseudoJetStructureBase>(s));
173 
174  return result_local;
175  } else {
176  // no subjets for seeding boosted tau reconstruction found, return an empty PseudoJet
177  if ( verbosity_ >= 1 ) {
178  std::cout << "No subjets found." << std::endl;
179  }
180  return PseudoJet();
181  }
182  }
void dumpSubJetStructure(const fastjet::PseudoJet &jet, int depth, int maxDepth, const std::string &depth_and_idx_string) const
int maxDepth_
the max depth for descending into clustering sequence
T sqrt(T t)
Definition: SSEVec.h:48
std::pair< PseudoJet, PseudoJet > findSubjets(const PseudoJet &jet, int depth, bool &subjetsFound) const
int verbosity_
flag to enable/disable debug output
static std::string join(char **cmd)
Definition: RemoteFile.cc:18
tuple cout
Definition: gather_cfg.py:121

Member Data Documentation

double contrib::CMSBoostedTauSeedingAlgorithm::dRMax_
private

the max value of the dR parameter

Definition at line 81 of file CMSBoostedTauSeedingAlgorithm.h.

double contrib::CMSBoostedTauSeedingAlgorithm::dRMin_
private

the min value of the dR parameter

Definition at line 80 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets().

int contrib::CMSBoostedTauSeedingAlgorithm::maxDepth_
private

the max depth for descending into clustering sequence

Definition at line 82 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets(), and result().

double contrib::CMSBoostedTauSeedingAlgorithm::muMax_
private

the max value of the mass-drop parameter

Definition at line 77 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets().

double contrib::CMSBoostedTauSeedingAlgorithm::muMin_
private

the min value of the mass-drop parameter

Definition at line 76 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets().

double contrib::CMSBoostedTauSeedingAlgorithm::ptMin_
private

minimum sub-jet pt

Definition at line 75 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets().

int contrib::CMSBoostedTauSeedingAlgorithm::verbosity_
private

flag to enable/disable debug output

Definition at line 84 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets(), and result().

double contrib::CMSBoostedTauSeedingAlgorithm::yMax_
private

the max value of the asymmetry parameter

Definition at line 79 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets().

double contrib::CMSBoostedTauSeedingAlgorithm::yMin_
private

the min value of the asymmetry parameter

Definition at line 78 of file CMSBoostedTauSeedingAlgorithm.h.

Referenced by findSubjets().