20 vector<fastjet::PseudoJet> & hardjetsOutput,
21 boost::shared_ptr<fastjet::ClusterSequence> & fjClusterSeq
24 if ( verbose_ )
cout <<
"Welcome to CATopSubJetAlgorithm::run" << endl;
31 for (
unsigned i = 0;
i < cell_particles.size(); ++
i) {
32 sumEt += cell_particles[
i].perp();
37 for (
unsigned int i = 0;
i < sumEtBins_.size(); ++
i ) {
38 if ( sumEt > sumEtBins_[
i] ) sumEtBinId =
i;
40 if ( verbose_ )
cout <<
"Using sumEt = " << sumEt <<
", bin = " << sumEtBinId << endl;
43 if ( sumEtBinId < 0 ) {
48 fastjet::PseudoJet blankJetA(0,0,0,0);
49 blankJetA.set_user_index(-1);
50 const fastjet::PseudoJet blankJet = blankJetA;
53 double deltarcut = deltarBins_[sumEtBinId];
54 double nCellMin = nCellBins_[sumEtBinId];
56 if ( verbose_ )
cout<<
"useAdjacency_ = "<<useAdjacency_<<endl;
57 if ( verbose_ && useAdjacency_==0)
cout<<
"No Adjacency"<<endl;
58 if ( verbose_ && useAdjacency_==1)
cout<<
"using deltar adjacency"<<endl;
59 if ( verbose_ && useAdjacency_==2)
cout<<
"using modified adjacency"<<endl;
60 if ( verbose_ && useAdjacency_==3)
cout<<
"using calorimeter nearest neighbor based adjacency"<<endl;
61 if ( verbose_ && useAdjacency_==1)
cout <<
"Using deltarcut = " << deltarcut << endl;
62 if ( verbose_ && useAdjacency_==3)
cout <<
"Using nCellMin = " << nCellMin << endl;
64 if ( verbose_ )
cout <<
"About to do jet clustering in CA" << endl;
78 if ( verbose_ )
cout <<
"Getting inclusive jets" << endl;
80 vector<fastjet::PseudoJet> inclusiveJets = fjClusterSeq->inclusive_jets(ptMin_);
82 if ( verbose_ )
cout <<
"Getting central jets" << endl;
84 vector<fastjet::PseudoJet> centralJets;
85 for (
unsigned int i = 0;
i < inclusiveJets.size();
i++) {
87 if (inclusiveJets[
i].
perp() > ptMin_ && fabs(inclusiveJets[
i].rapidity()) < centralEtaCut_) {
88 centralJets.push_back(inclusiveJets[
i]);
93 sort( centralJets.begin(), centralJets.end(), compEt );
96 vector<math::XYZTLorentzVector> p4_hardJets;
100 vector<vector<int> > indices( centralJets.size() );
103 vector<fastjet::PseudoJet>::iterator jetIt = centralJets.begin(),
104 centralJetsEnd = centralJets.end();
105 if ( verbose_ )
cout<<
"Loop over jets"<<endl;
107 for ( ; jetIt != centralJetsEnd; ++jetIt ) {
108 if ( verbose_ )
cout<<
"\nJet "<<i<<endl;
110 fastjet::PseudoJet localJet = *jetIt;
116 if ( verbose_ )
cout<<
"local jet pt = "<<localJet.perp()<<endl;
117 if ( verbose_ )
cout<<
"deltap = "<<ptFracBins_[sumEtBinId]<<endl;
119 double ptHard = ptFracBins_[sumEtBinId]*localJet.perp();
120 vector<fastjet::PseudoJet> leftoversAll;
123 if ( verbose_ )
cout <<
"Doing decomposition 1" << endl;
124 fastjet::PseudoJet ja, jb;
125 vector<fastjet::PseudoJet> leftovers1;
126 bool hardBreak1 = decomposeJet(localJet,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,ja,jb,leftovers1);
127 leftoversAll.insert(leftoversAll.end(),leftovers1.begin(),leftovers1.end());
132 if ( verbose_ )
cout <<
"Doing decomposition 2. ja->jaa+jab?" << endl;
133 fastjet::PseudoJet jaa, jab;
134 vector<fastjet::PseudoJet> leftovers2a;
135 bool hardBreak2a =
false;
136 if (hardBreak1) hardBreak2a = decomposeJet(ja,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,jaa,jab,leftovers2a);
137 leftoversAll.insert(leftoversAll.end(),leftovers2a.begin(),leftovers2a.end());
139 if ( verbose_ )
cout <<
"Doing decomposition 2. ja->jba+jbb?" << endl;
140 fastjet::PseudoJet jba, jbb;
141 vector<fastjet::PseudoJet> leftovers2b;
142 bool hardBreak2b =
false;
143 if (hardBreak1) hardBreak2b = decomposeJet(jb,*fjClusterSeq,cell_particles,ptHard,nCellMin,deltarcut,jba,jbb,leftovers2b);
144 leftoversAll.insert(leftoversAll.end(),leftovers2b.begin(),leftovers2b.end());
151 if ( verbose_ )
cout <<
"Done with decomposition" << endl;
153 if ( verbose_ )
cout<<
"hardBreak1 = "<<hardBreak1<<endl;
154 if ( verbose_ )
cout<<
"hardBreak2a = "<<hardBreak2a<<endl;
155 if ( verbose_ )
cout<<
"hardBreak2b = "<<hardBreak2b<<endl;
157 fastjet::PseudoJet hardA = blankJet, hardB = blankJet, hardC = blankJet, hardD = blankJet;
163 if(verbose_)
cout<<
"Hardbreak failed. Save subjet1=localJet"<<endl;
165 if (hardBreak1 && !hardBreak2a && !hardBreak2b) {
170 if(verbose_)
cout<<
"First decomposition succeeded, both second decompositions failed. Save subjet1=ja subjet2=jb"<<endl;
172 if (hardBreak1 && hardBreak2a && !hardBreak2b) {
177 if(verbose_)
cout<<
"First decomposition succeeded, ja split succesfully, jb did not split. Save subjet1=jaa subjet2=jab subjet3=jb"<<endl;
179 if (hardBreak1 && !hardBreak2a && hardBreak2b) {
184 if(verbose_)
cout<<
"First decomposition succeeded, jb split succesfully, ja did not split. Save subjet1=jba subjet2=jbb subjet3=ja"<<endl;
186 if (hardBreak1 && hardBreak2a && hardBreak2b) {
191 if(verbose_)
cout<<
"First decomposition and both secondary decompositions succeeded. Save subjet1=jaa subjet2=jab subjet3=jba subjet4=jbb"<<endl;
195 fastjet::PseudoJet subjet1 = blankJet;
196 fastjet::PseudoJet subjet2 = blankJet;
197 fastjet::PseudoJet subjet3 = blankJet;
198 fastjet::PseudoJet subjet4 = blankJet;
199 subjet1 = hardA; subjet2 = hardB; subjet3 = hardC; subjet4 = hardD;
202 vector<fastjet::PseudoJet> hardSubjets;
205 std::cout <<
"HardA : user_index = " << hardA.user_index() <<
", (Pt,Y,Phi,M) = ("
206 << hardA.pt() <<
", " << hardA.rapidity() <<
", "
207 << hardA.phi() <<
", " << hardA.m() <<
")" << std::endl;
209 std::cout <<
"HardB : user_index = " << hardB.user_index() <<
", (Pt,Y,Phi,M) = ("
210 << hardB.pt() <<
", " << hardB.rapidity() <<
", "
211 << hardB.phi() <<
", " << hardB.m() <<
")" << std::endl;
213 std::cout <<
"HardC : user_index = " << hardC.user_index() <<
", (Pt,Y,Phi,M) = ("
214 << hardC.pt() <<
", " << hardC.rapidity() <<
", "
215 << hardC.phi() <<
", " << hardC.m() <<
")" << std::endl;
217 std::cout <<
"HardD : user_index = " << hardD.user_index() <<
", (Pt,Y,Phi,M) = ("
218 << hardD.pt() <<
", " << hardD.rapidity() <<
", "
219 << hardD.phi() <<
", " << hardD.m() <<
")" << std::endl;
226 if ( subjet1.pt() > 0.0001 )
227 hardSubjets.push_back(subjet1);
228 if ( subjet2.pt() > 0.0001 )
229 hardSubjets.push_back(subjet2);
230 if ( subjet3.pt() > 0.0001 )
231 hardSubjets.push_back(subjet3);
232 if ( subjet4.pt() > 0.0001 )
233 hardSubjets.push_back(subjet4);
234 sort(hardSubjets.begin(), hardSubjets.end(), compEt );
237 fastjet::PseudoJet candidate =
join(hardSubjets);
239 candidate.reset_momentum( jetIt->px(), jetIt->py(), jetIt->pz(), jetIt->e() );
242 std::cout <<
"Final top-jet candidate: (Pt,Y,Phi,M) = ("
243 << candidate.pt() <<
", " << candidate.rapidity() <<
", "
244 << candidate.phi() <<
", " << candidate.m() <<
")" << std::endl;
245 std::vector<fastjet::PseudoJet>
pieces = candidate.pieces();
246 std::cout <<
"Number of pieces = " << pieces.size() << std::endl;
247 for ( std::vector<fastjet::PseudoJet>::const_iterator ibegin = pieces.begin(), iend = pieces.end(), i = ibegin;
249 std::cout <<
" Piece : " << i - ibegin <<
", (Pt,Y,Phi,M) = ("
250 << i->pt() <<
", " << i->rapidity() <<
", "
251 << i->phi() <<
", " << i->m() <<
")" << std::endl;
255 hardjetsOutput.push_back( candidate );
269 const vector<fastjet::PseudoJet> & cell_particles,
270 const fastjet::ClusterSequence & theClusterSequence,
271 double nCellMin )
const {
274 double eta1 = jet1.rapidity();
275 double phi1 = jet1.phi();
276 double eta2 = jet2.rapidity();
277 double phi2 = jet2.phi();
279 double deta =
abs(eta2 - eta1) / 0.087;
282 return ( ( deta + dphi ) <= nCellMin );
291 const fastjet::ClusterSequence & theClusterSequence,
292 const vector<fastjet::PseudoJet> & cell_particles,
293 double ptHard,
double nCellMin,
double deltarcut,
294 fastjet::PseudoJet & ja, fastjet::PseudoJet & jb,
295 vector<fastjet::PseudoJet> & leftovers)
const {
298 fastjet::PseudoJet
j = theJet;
299 double InputObjectPt = j.perp();
300 if ( verbose_ )
cout<<
"Input Object Pt = "<<InputObjectPt<<endl;
301 if ( verbose_ )
cout<<
"ptHard = "<<ptHard<<endl;
303 if ( verbose_ )
cout<<
"start while loop"<<endl;
306 goodBreak = theClusterSequence.has_parents(j,ja,jb);
308 if ( verbose_ )
cout<<
"bad break. this is one cell. can't decluster anymore."<<endl;
312 if ( verbose_ )
cout<<
"good break. ja Pt = "<<ja.perp()<<
" jb Pt = "<<jb.perp()<<endl;
317 double clusters_deltar=fabs(ja.eta()-jb.eta())+fabs(
deltaPhi(ja.phi(),jb.phi()));
319 if ( verbose_ && useAdjacency_ ==1)
cout<<
"clusters_deltar = "<<clusters_deltar<<endl;
320 if ( verbose_ && useAdjacency_ ==1)
cout<<
"deltar cut = "<<deltarcut<<endl;
322 if ( useAdjacency_==1 && clusters_deltar < deltarcut){
323 if ( verbose_ )
cout<<
"clusters too close. consant adj. break."<<endl;
328 double clusters_deltaR=
deltaR( ja.rapidity(), ja.phi(), jb.rapidity(), jb.phi() );
330 if ( verbose_ && useAdjacency_ ==2)
cout<<
"clusters_deltaR = "<<clusters_deltaR<<endl;
331 if ( verbose_ && useAdjacency_ ==2)
cout<<
"0.4-0.0004*InputObjectPt = "<<0.4-0.0004*InputObjectPt<<endl;
333 if ( useAdjacency_==2 && clusters_deltaR < 0.4-0.0004*InputObjectPt)
335 if ( verbose_ )
cout<<
"clusters too close. modified adj. break."<<endl;
340 if ( useAdjacency_==3 && adjacentCells(ja,jb,cell_particles,theClusterSequence,nCellMin) ){
341 if ( verbose_ )
cout<<
"clusters too close in the calorimeter. calorimeter adj. break."<<endl;
345 if ( verbose_ )
cout<<
"clusters pass distance cut"<<endl;
349 if ( verbose_ )
cout<<
"ptHard = "<<ptHard<<endl;
351 if (ja.perp() < ptHard && jb.perp() < ptHard){
352 if ( verbose_ )
cout<<
"two soft clusters. dead end"<<endl;
356 if (ja.perp() > ptHard && jb.perp() > ptHard){
357 if ( verbose_ )
cout<<
"two hard clusters. done"<<endl;
361 else if (ja.perp() > jb.perp()) {
362 if ( verbose_ )
cout<<
"ja hard jb soft. try to split hard. j = ja"<<endl;
364 vector<fastjet::PseudoJet> particles = theClusterSequence.constituents(jb);
365 leftovers.insert(leftovers.end(),particles.begin(),particles.end());
368 if ( verbose_ )
cout<<
"ja hard jb soft. try to split hard. j = jb"<<endl;
370 vector<fastjet::PseudoJet> particles = theClusterSequence.constituents(ja);
371 leftovers.insert(leftovers.end(),particles.begin(),particles.end());
375 if ( verbose_ )
cout<<
"did not decluster."<<endl;
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
bool adjacentCells(const fastjet::PseudoJet &jet1, const fastjet::PseudoJet &jet2, const std::vector< fastjet::PseudoJet > &cell_particles, const fastjet::ClusterSequence &theClusterSequence, double nCellMin) const
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
Abs< T >::type abs(const T &t)
double deltaR(double eta1, double eta2, double phi1, double phi2)
double deltaPhi(double phi1, double phi2)
static std::string join(char **cmd)
T perp() const
Magnitude of transverse component.
void run(const std::vector< fastjet::PseudoJet > &cell_particles, std::vector< fastjet::PseudoJet > &hardjetsOutput, boost::shared_ptr< fastjet::ClusterSequence > &fjClusterSeq)
Find the ProtoJets from the collection of input Candidates.