21 #include <fastjet/ClusterSequenceArea.hh> 32 ostream &
operator<<(ostream & ostr, fastjet::PseudoJet &
jet);
51 int activeAreaRepeats,
54 : moduleLabel_(moduleLabel)
55 , jetAlgorithm_(jetAlgorithm)
60 , massDropCut_(massDropCut)
61 , asymmCut2_(asymmCut*asymmCut)
62 , asymmCutLater_(asymmCutLater)
63 , doAreaFastjet_(doAreaFastjet)
64 , ghostEtaMax_(ghostEtaMax)
65 , activeAreaRepeats_(activeAreaRepeats)
66 , ghostArea_(ghostArea)
83 <<
"Jet Algorithm for SubjetFilterAlgorithm is invalid: " 84 <<
jetAlgorithm_<<
", use (ca|CambridgeAachen)|(Kt|kt)|(AntiKt|ak)"<<endl;
89 fjAreaDef_=
new fastjet::AreaDefinition(fastjet::active_area_explicit_ghosts,
109 std::vector<CompoundPseudoJet>& fjJets,
118 new fastjet::ClusterSequence (fjInputs,*
fjJetDef_);
120 vector<fastjet::PseudoJet> fjFatJets =
121 fastjet::sorted_by_pt(cs->inclusive_jets(
jetPtMin_));
126 for (
size_t iFat=0;iFat<nFat;iFat++) {
128 if (
verbose_)
cout<<endl<<iFat<<
". FATJET: "<<fjFatJets[iFat]<<endl;
130 fastjet::PseudoJet fjFatJet = fjFatJets[iFat];
131 fastjet::PseudoJet fjCurrentJet(fjFatJet);
132 fastjet::PseudoJet fjSubJet1,fjSubJet2;
135 vector<CompoundPseudoSubJet> subJets;
139 while ((hadSubJets = cs->has_parents(fjCurrentJet,fjSubJet1,fjSubJet2))) {
141 if (fjSubJet1.m() < fjSubJet2.m())
swap(fjSubJet1,fjSubJet2);
143 if (
verbose_)
cout<<
"SUBJET CANDIDATES: "<<fjSubJet1<<endl
144 <<
" "<<fjSubJet2<<endl;
145 if (
verbose_)
cout<<
"md="<<fjSubJet1.m()/fjCurrentJet.m()
146 <<
" y="<<fjSubJet1.kt_distance(fjSubJet2)/fjCurrentJet.m2()
151 fjSubJet1.kt_distance(fjSubJet2)>
asymmCut2_*fjCurrentJet.m2())) {
155 fjCurrentJet = fjSubJet1;
167 fjSubJet1.kt_distance(fjSubJet2)<=
asymmCut2_*fjCurrentJet.m2()) {
174 vector<fastjet::PseudoJet> fjFilterJets;
175 double Rbb =
std::sqrt(fjSubJet1.squared_distance(fjSubJet2));
178 fjFilterJets=fastjet::sorted_by_pt(cs->exclusive_subjets(fjCurrentJet,dcut));
181 cout<<
"Rbb="<<Rbb<<
", Rfilt="<<Rfilt<<endl;
182 cout<<
"FILTER JETS: "<<flush;
183 for (
size_t i=0;
i<fjFilterJets.size();
i++) {
184 if (
i>0)
cout<<
" "<<flush;
185 cout<<fjFilterJets[
i]<<endl;
189 vector<fastjet::PseudoJet> fjSubJets;
190 fjSubJets.push_back(fjSubJet1);
191 fjSubJets.push_back(fjSubJet2);
192 size_t nFilter = fjFilterJets.size();
193 for (
size_t iFilter=0;iFilter<nFilter;iFilter++)
194 fjSubJets.push_back(fjFilterJets[iFilter]);
196 for (
size_t iSub=0;iSub<fjSubJets.size();iSub++) {
197 vector<fastjet::PseudoJet> fjConstituents=
198 cs->constituents(fjSubJets[iSub]);
199 vector<int> constituents;
200 for (
size_t iConst=0;iConst<fjConstituents.size();iConst++) {
201 int userIndex = fjConstituents[iConst].user_index();
202 if (userIndex>=0) constituents.push_back(userIndex);
206 ((fastjet::ClusterSequenceArea*)
cs)->
area(fjSubJets[iSub]) : 0.0;
208 if (iSub<2||!constituents.empty())
218 if (
verbose_)
cout<<
"write fatjet with "<<subJets.size()<<
" sub+filter jets" 222 ((fastjet::ClusterSequenceArea*)
cs)->
area(fjFatJet) : 0.0;
227 if (subJets.size()>3)
nfound_++;
231 if (
verbose_)
cout<<endl<<fjJets.size()<<
" FATJETS written\n"<<endl;
243 std::stringstream ss;
244 ss<<
"************************************************************\n" 245 <<
"* "<<
moduleLabel_<<
" (SubjetFilterAlgorithm) SUMMARY:\n" 246 <<
"************************************************************\n" 250 <<
"eff = "<<eff<<endl
251 <<
"************************************************************\n";
259 ostr <<
"pt=" <<setw(10)<<jet.perp()
260 <<
" eta="<<setw(6) <<jet.eta()
261 <<
" m=" <<setw(10)<<jet.m();
unique_ptr< ClusterSequence > cs
fastjet::AreaDefinition * fjAreaDef_
void swap(Association< C > &lhs, Association< C > &rhs)
void run(const std::vector< fastjet::PseudoJet > &inputs, std::vector< CompoundPseudoJet > &fatJets, const edm::EventSetup &iSetup)
CompoundPseudoJet holds an association of fastjet::PseudoJets that represent a "hard" top jet with su...
fastjet::JetDefinition * fjJetDef_
SubjetFilterAlgorithm(const std::string &moduleLabel, const std::string &jetAlgorithm, unsigned nFatMax, double rParam, double rFilt, double jetPtMin, double massDropCut, double asymmCut, bool asymmCutLater, bool doAreaFastjet, double ghostEtaMax, int activeAreaRepeats, double ghostArea, bool verbose)
ostream & operator<<(ostream &ostr, fastjet::PseudoJet &jet)
does the actual work for printing out a jet
std::ostream & operator<<(std::ostream &out, const std::tuple< Types... > &value)
virtual ~SubjetFilterAlgorithm()
std::string jetAlgorithm_
std::string summary() const