Find the ProtoJets from the collection of input Candidates.
25 if (
verbose_ )
cout <<
"Welcome to CATopSubJetAlgorithm::run" << endl;
32 for (
unsigned i = 0;
i < cell_particles.size(); ++
i) {
33 sumEt += cell_particles[
i].perp();
41 if (
verbose_ )
cout <<
"Using sumEt = " << sumEt <<
", bin = " << sumEtBinId << endl;
44 if ( sumEtBinId < 0 ) {
49 fastjet::PseudoJet blankJetA(0,0,0,0);
50 blankJetA.set_user_index(-1);
51 const fastjet::PseudoJet blankJet = blankJetA;
65 if (
verbose_ )
cout <<
"About to do jet clustering in CA" << endl;
79 if (
verbose_ )
cout <<
"Getting inclusive jets" << endl;
81 vector<fastjet::PseudoJet> inclusiveJets = fjClusterSeq->inclusive_jets(
ptMin_);
85 vector<fastjet::PseudoJet> centralJets;
86 for (
unsigned int i = 0;
i < inclusiveJets.size();
i++) {
89 centralJets.push_back(inclusiveJets[
i]);
94 sort( centralJets.begin(), centralJets.end(), compEt );
97 vector<math::XYZTLorentzVector> p4_hardJets;
101 vector<vector<int> >
indices( centralJets.size() );
104 vector<fastjet::PseudoJet>::iterator jetIt = centralJets.begin(),
105 centralJetsEnd = centralJets.end();
108 for ( ; jetIt != centralJetsEnd; ++jetIt ) {
111 fastjet::PseudoJet localJet = *jetIt;
117 if (
verbose_ )
cout<<
"local jet pt = "<<localJet.perp()<<endl;
120 double ptHard =
ptFracBins_[sumEtBinId]*localJet.perp();
121 vector<fastjet::PseudoJet> leftoversAll;
124 if (
verbose_ )
cout <<
"Doing decomposition 1" << endl;
125 fastjet::PseudoJet ja, jb;
126 vector<fastjet::PseudoJet> leftovers1;
127 bool hardBreak1 =
decomposeJet(localJet,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,ja,jb,leftovers1);
128 leftoversAll.insert(leftoversAll.end(),leftovers1.begin(),leftovers1.end());
133 if (
verbose_ )
cout <<
"Doing decomposition 2. ja->jaa+jab?" << endl;
134 fastjet::PseudoJet jaa, jab;
135 vector<fastjet::PseudoJet> leftovers2a;
136 bool hardBreak2a =
false;
137 if (hardBreak1) hardBreak2a =
decomposeJet(ja,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,jaa,jab,leftovers2a);
138 leftoversAll.insert(leftoversAll.end(),leftovers2a.begin(),leftovers2a.end());
140 if (
verbose_ )
cout <<
"Doing decomposition 2. ja->jba+jbb?" << endl;
141 fastjet::PseudoJet jba, jbb;
142 vector<fastjet::PseudoJet> leftovers2b;
143 bool hardBreak2b =
false;
144 if (hardBreak1) hardBreak2b =
decomposeJet(jb,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,jba,jbb,leftovers2b);
145 leftoversAll.insert(leftoversAll.end(),leftovers2b.begin(),leftovers2b.end());
152 if (
verbose_ )
cout <<
"Done with decomposition" << endl;
154 if (
verbose_ )
cout<<
"hardBreak1 = "<<hardBreak1<<endl;
155 if (
verbose_ )
cout<<
"hardBreak2a = "<<hardBreak2a<<endl;
156 if (
verbose_ )
cout<<
"hardBreak2b = "<<hardBreak2b<<endl;
158 fastjet::PseudoJet hardA = blankJet, hardB = blankJet, hardC = blankJet, hardD = blankJet;
164 if(
verbose_)
cout<<
"Hardbreak failed. Save subjet1=localJet"<<endl;
166 if (hardBreak1 && !hardBreak2a && !hardBreak2b) {
171 if(
verbose_)
cout<<
"First decomposition succeeded, both second decompositions failed. Save subjet1=ja subjet2=jb"<<endl;
173 if (hardBreak1 && hardBreak2a && !hardBreak2b) {
178 if(
verbose_)
cout<<
"First decomposition succeeded, ja split succesfully, jb did not split. Save subjet1=jaa subjet2=jab subjet3=jb"<<endl;
180 if (hardBreak1 && !hardBreak2a && hardBreak2b) {
185 if(
verbose_)
cout<<
"First decomposition succeeded, jb split succesfully, ja did not split. Save subjet1=jba subjet2=jbb subjet3=ja"<<endl;
187 if (hardBreak1 && hardBreak2a && hardBreak2b) {
192 if(
verbose_)
cout<<
"First decomposition and both secondary decompositions succeeded. Save subjet1=jaa subjet2=jab subjet3=jba subjet4=jbb"<<endl;
196 fastjet::PseudoJet subjet1 = blankJet;
197 fastjet::PseudoJet subjet2 = blankJet;
198 fastjet::PseudoJet subjet3 = blankJet;
199 fastjet::PseudoJet subjet4 = blankJet;
200 subjet1 = hardA; subjet2 = hardB; subjet3 = hardC; subjet4 = hardD;
203 vector<fastjet::PseudoJet> hardSubjets;
206 std::cout <<
"HardA : user_index = " << hardA.user_index() <<
", (Pt,Y,Phi,M) = ("
207 << hardA.pt() <<
", " << hardA.rapidity() <<
", "
208 << hardA.phi() <<
", " << hardA.m() <<
")" << std::endl;
210 std::cout <<
"HardB : user_index = " << hardB.user_index() <<
", (Pt,Y,Phi,M) = ("
211 << hardB.pt() <<
", " << hardB.rapidity() <<
", "
212 << hardB.phi() <<
", " << hardB.m() <<
")" << std::endl;
214 std::cout <<
"HardC : user_index = " << hardC.user_index() <<
", (Pt,Y,Phi,M) = ("
215 << hardC.pt() <<
", " << hardC.rapidity() <<
", "
216 << hardC.phi() <<
", " << hardC.m() <<
")" << std::endl;
218 std::cout <<
"HardD : user_index = " << hardD.user_index() <<
", (Pt,Y,Phi,M) = ("
219 << hardD.pt() <<
", " << hardD.rapidity() <<
", "
220 << hardD.phi() <<
", " << hardD.m() <<
")" << std::endl;
227 if ( subjet1.pt() > 0.0001 )
228 hardSubjets.push_back(subjet1);
229 if ( subjet2.pt() > 0.0001 )
230 hardSubjets.push_back(subjet2);
231 if ( subjet3.pt() > 0.0001 )
232 hardSubjets.push_back(subjet3);
233 if ( subjet4.pt() > 0.0001 )
234 hardSubjets.push_back(subjet4);
235 sort(hardSubjets.begin(), hardSubjets.end(), compEt );
238 fastjet::PseudoJet candidate =
join(hardSubjets);
240 candidate.reset_momentum( jetIt->px(), jetIt->py(), jetIt->pz(), jetIt->e() );
243 std::cout <<
"Final top-jet candidate: (Pt,Y,Phi,M) = ("
244 << candidate.pt() <<
", " << candidate.rapidity() <<
", "
245 << candidate.phi() <<
", " << candidate.m() <<
")" << std::endl;
246 std::vector<fastjet::PseudoJet>
pieces = candidate.pieces();
247 std::cout <<
"Number of pieces = " << pieces.size() << std::endl;
248 for ( std::vector<fastjet::PseudoJet>::const_iterator ibegin = pieces.begin(), iend = pieces.end(), i = ibegin;
250 std::cout <<
" Piece : " << i - ibegin <<
", (Pt,Y,Phi,M) = ("
251 << i->pt() <<
", " << i->rapidity() <<
", "
252 << i->phi() <<
", " << i->m() <<
")" << std::endl;
256 hardjetsOutput.push_back( candidate );
std::vector< double > sumEtBins_
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
bool decomposeJet(const fastjet::PseudoJet &theJet, const fastjet::ClusterSequence &theClusterSequence, const std::vector< fastjet::PseudoJet > &cell_particles, double ptHard, double nCellMin, double deltarcut, fastjet::PseudoJet &ja, fastjet::PseudoJet &jb, std::vector< fastjet::PseudoJet > &leftovers) const
std::vector< double > ptFracBins_
std::vector< double > deltarBins_
static std::string join(char **cmd)
T perp() const
Magnitude of transverse component.
std::vector< double > nCellBins_