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 | Private Member Functions | Private Attributes
CMSMidpointAlgorithm Class Reference

#include <CMSMidpointAlgorithm.h>

Public Types

typedef std::vector< ProtoJet * > InternalCollection
 

Public Member Functions

 CMSMidpointAlgorithm ()
 Default constructor. More...
 
 CMSMidpointAlgorithm (double fSeedThreshold, double fConeRadius, double fConeAreaFraction, int fMaxPairSize, int fMaxIterations, double fOverlapThreshold, int fDebugLevel)
 
void run (const JetReco::InputCollection &fInput, JetReco::OutputCollection *fOutput)
 

Private Member Functions

void addClustersToPairs (const JetReco::InputCollection &fInput, std::vector< int > &testPair, std::vector< std::vector< int > > &pairs, std::vector< std::vector< bool > > &distanceOK, int maxClustersInPair)
 
void findStableConesFromMidPoints (const JetReco::InputCollection &fInput, InternalCollection *fOutput)
 
void findStableConesFromSeeds (const JetReco::InputCollection &fInput, InternalCollection *fOutput)
 
void iterateCone (const JetReco::InputCollection &fInput, double startRapidity, double startPhi, double startPt, bool reduceConeSize, InternalCollection *fOutput)
 
void splitAndMerge (const JetReco::InputCollection &fInput, InternalCollection *fProtoJets, JetReco::OutputCollection *fFinalJets)
 

Private Attributes

double theConeAreaFraction
 
double theConeRadius
 
int theDebugLevel
 
int theMaxIterations
 
int theMaxPairSize
 
double theOverlapThreshold
 
double theSeedThreshold
 

Detailed Description

CMSMidpointAlgorithm is an algorithm for CMS jet reconstruction baded on the CDF midpoint algorithm.

The algorithm is documented in the proceedings of the Physics at RUN II: QCD and Weak Boson Physics Workshop: hep-ex/0005012.

The algorithm runs off of generic Candidates

Author
Robert M Harris, Fermilab
Version
1st Version Feb. 4, 2005 Based on the CDF Midpoint Algorithm code by Matthias Toennesmann.
2nd Version Apr. 6, 2005 Modifications toward integration in new EDM.
3rd Version Oct. 19, 2005 Modified to work with real CaloTowers from Jeremy Mans
F.Ratnikov, Mar. 8, 2006. Work from Candidate

Definition at line 26 of file CMSMidpointAlgorithm.h.

Member Typedef Documentation

Definition at line 29 of file CMSMidpointAlgorithm.h.

Constructor & Destructor Documentation

CMSMidpointAlgorithm::CMSMidpointAlgorithm ( )
inline

Default constructor.

Definition at line 32 of file CMSMidpointAlgorithm.h.

CMSMidpointAlgorithm::CMSMidpointAlgorithm ( double  fSeedThreshold,
double  fConeRadius,
double  fConeAreaFraction,
int  fMaxPairSize,
int  fMaxIterations,
double  fOverlapThreshold,
int  fDebugLevel 
)
inline
Constructor takes as input all the values of the algorithm that the user can change.
Parameters
fSeedThreshold,:minimum ET in GeV of an CaloTower that can seed a jet.
fTowerThreshold,:minimum ET in GeV of an CaloTower that is included in a jet
fConeRadius,:nominal radius of the jet in eta-phi space
fConeAreaFraction,:multiplier to reduce the search cone area during the iteration phase. introduced by CDF in 2002 to avoid energy loss due to proto-jet migration. Actively being discussed in 2005. A value of 1.0 gives the original run II algorithm in hep-ex/0005012. CDF value of 0.25 gives new search cone of 0.5 theConeRadius during iteration, but keeps the final cone at theConeRadius for the last iteration.
fMaxPairSize,:Maximum size of proto-jet pair, triplet, etc for defining midpoint. Both CDF and D0 use 2.
fMaxIterations,:Maximum number of iterations before finding a stable cone.
fOverlapThreshold,:When two proto-jets overlap, this is the merging threshold on the fraction of PT in the overlap region compared to the lower Pt Jet. If the overlapPt/lowerJetPt > theOverlapThreshold the 2 proto-jets will be merged into one final jet. if the overlapPt/lowerJetPt < theOverlapThreshold the towers in the two proto-jets will be seperated into two distinct sets of towers: two final jets. D0 has used 0.5, CDF has used both 0.5 and 0.75 in run 2, and 0.75 in run 1.
fDebugLevel,:integer level of diagnostic printout: 0 = no printout, 1 = minimal printout, 2 = pages per event.

Definition at line 64 of file CMSMidpointAlgorithm.h.

65  :
66  theSeedThreshold(fSeedThreshold),
67  theConeRadius(fConeRadius),
68  theConeAreaFraction(fConeAreaFraction),
69  theMaxPairSize(fMaxPairSize),
70  theMaxIterations(fMaxIterations),
71  theOverlapThreshold(fOverlapThreshold),
72  theDebugLevel(fDebugLevel)
73  { }

Member Function Documentation

void CMSMidpointAlgorithm::addClustersToPairs ( const JetReco::InputCollection fInput,
std::vector< int > &  testPair,
std::vector< std::vector< int > > &  pairs,
std::vector< std::vector< bool > > &  distanceOK,
int  maxClustersInPair 
)
private

Add proto-jets to pairs, triplets, etc, prior to finding their midpoints. Called by findStableConesFromMidPoints but public method to allow testing and studies.

Definition at line 278 of file CMSMidpointAlgorithm.cc.

281 {
282  // Recursively adds clusters to pairs, triplets, ... whose mid-points are then calculated.
283 
284  // Find StableCone number to start with (either 0 at the beginning or last element of testPair + 1).
285  int nextClusterStart = 0;
286  if(testPair.size())
287  nextClusterStart = testPair.back() + 1;
288  for(unsigned int nextCluster = nextClusterStart; nextCluster <= distanceOK.size(); ++nextCluster){
289  // Is new SeedCone less than 2*_theConeRadius apart from all clusters in testPair?
290  bool addCluster = true;
291  for(unsigned int iCluster = 0; iCluster < testPair.size() && addCluster; ++iCluster)
292  if(!distanceOK[nextCluster - 1][testPair[iCluster]])
293  addCluster = false;
294  if(addCluster){
295  // Add it to the testPair.
296  testPair.push_back(nextCluster);
297  // If testPair is a pair, add it to pairs.
298  if(testPair.size() > 1)
299  pairs.push_back(testPair);
300  // If not bigger than allowed, find more clusters within 2*theConeRadius.
301  if(testPair.size() < unsigned(maxClustersInPair))
302  addClustersToPairs(fInput, testPair,pairs,distanceOK,maxClustersInPair);
303  // All combinations containing testPair found. Remove last element.
304  testPair.pop_back();
305  }
306  }
307 }
void addClustersToPairs(const JetReco::InputCollection &fInput, std::vector< int > &testPair, std::vector< std::vector< int > > &pairs, std::vector< std::vector< bool > > &distanceOK, int maxClustersInPair)
void CMSMidpointAlgorithm::findStableConesFromMidPoints ( const JetReco::InputCollection fInput,
InternalCollection fOutput 
)
private

Add to the list of proto-jets the list of midpoints. Called by run but public method to allow testing and studies.

Definition at line 236 of file CMSMidpointAlgorithm.cc.

References gather_cfg::cout, Geom::deltaR2(), ProtoJet::phi(), python.multivaluedict::sort(), and ProtoJet::y().

236  {
237  // We take the previous list of stable protojets from seeds as input and add to it
238  // those from the midpoints between proto-jet pairs, triplets, etc.
239  // distanceOK[i-1][j] = Is distance between stableCones i and j (i>j) less than 2*_theConeRadius?
240  vector< vector<bool> > distanceOK; // A vector of vectors
241  distanceOK.resize(stableCones->size() - 1); // Set the outer vector size, num protojets - 1
242  for(unsigned int nCluster1 = 1; nCluster1 < stableCones->size(); ++nCluster1){ // Loop over the protojets
243  distanceOK[nCluster1 - 1].resize(nCluster1); //Set inner vector size: monotonically increasing.
244  const ProtoJet* cluster1 = (*stableCones)[nCluster1];
245  for(unsigned int nCluster2 = 0; nCluster2 < nCluster1; ++nCluster2){ // Loop over the other proto-jets
246  const ProtoJet* cluster2 = (*stableCones)[nCluster2];
247  double dR2 = deltaR2 (cluster1->y(), cluster1->phi(), cluster2->y(), cluster2->phi());
248  distanceOK[nCluster1 - 1][nCluster2] = dR2 < 4*theConeRadius*theConeRadius;
249  }
250  }
251 
252  // Find all pairs (triplets, ...) of stableCones which are less than 2*theConeRadius apart from each other.
253  vector< vector<int> > pairs(0);
254  vector<int> testPair(0);
255  int maxClustersInPair = theMaxPairSize; // Set maximum proto-jets to a pair or a triplet, etc.
256  if(!maxClustersInPair)maxClustersInPair = stableCones->size(); // If zero, then skys the limit!
257  addClustersToPairs(fInput, testPair,pairs,distanceOK,maxClustersInPair); // Make the pairs, triplets, etc.
258 
259  // Loop over all combinations. Calculate MidPoint. Make midPointClusters.
260  bool reduceConeSize = false; // Note that here we keep the iteration cone size fixed.
261  for(vector<vector<int> >::const_iterator iPair = pairs.begin(); iPair != pairs.end(); ++iPair) {
262  const vector<int> & Pair = *iPair;
263  // Calculate rapidity, phi and energy of MidPoint.
264  reco::Particle::LorentzVector midPoint (0,0,0,0);
265  for(vector<int>::const_iterator iPairMember = Pair.begin(); iPairMember != Pair.end(); ++iPairMember) {
266  midPoint += (*stableCones)[*iPairMember]->p4();
267  }
268  if(theDebugLevel>=2)cout << endl << "[CMSMidpointAlgorithm] midpoint " << iPair-pairs.begin() << ": y = " << midPoint.Rapidity() << ", phi=" << midPoint.Phi() <<
269  ", size=" << Pair.size() << endl;
270  iterateCone(fInput, midPoint.Rapidity(),midPoint.Phi(),midPoint.e(),reduceConeSize,stableCones);
271  }
272  GreaterByPtPtr<ProtoJet> compJets;
273  sort (stableCones->begin(), stableCones->end(), compJets);
274 }
void iterateCone(const JetReco::InputCollection &fInput, double startRapidity, double startPhi, double startPt, bool reduceConeSize, InternalCollection *fOutput)
Transient Jet class used by the reconstruction algorithms.
Definition: ProtoJet.h:25
double y() const
Definition: ProtoJet.h:85
double phi() const
Definition: ProtoJet.h:81
double deltaR2(const Vector1 &v1, const Vector2 &v2)
Definition: VectorUtil.h:78
tuple cout
Definition: gather_cfg.py:121
void addClustersToPairs(const JetReco::InputCollection &fInput, std::vector< int > &testPair, std::vector< std::vector< int > > &pairs, std::vector< std::vector< bool > > &distanceOK, int maxClustersInPair)
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:27
void CMSMidpointAlgorithm::findStableConesFromSeeds ( const JetReco::InputCollection fInput,
InternalCollection fOutput 
)
private

Find the list of proto-jets from the seeds. Called by run, but public method to allow testing and studies.

Definition at line 130 of file CMSMidpointAlgorithm.cc.

References gather_cfg::cout, and i.

130  {
131  // This dictates that the cone size will be reduced in the iterations procedure,
132  // to prevent excessive cone movement, and then will be enlarged at the end of
133  // the iteration procedure (ala CDF).
134  bool reduceConeSize = true;
135 
136  // Get the Seed Towers sorted by Et.
137  vector<InputItem> seedTowers = etOrderedCaloTowers(fInput, theSeedThreshold); //This gets towers
138 
139  // Loop over all Seeds
140  for(vector<InputItem>::const_iterator i = seedTowers.begin(); i != seedTowers.end(); ++i) {
141  double seedEta = (*i)->eta ();
142  double seedPhi = (*i)->phi ();
143 
144  // Find stable cone from seed.
145  // This iterates the cone centroid, makes a proto-jet, and adds it to the list.
146  if(theDebugLevel>=2) cout << endl << "[CMSMidpointAlgorithm] seed " << i-seedTowers.begin() << ": eta=" << seedEta << ", phi=" << seedPhi << endl;
147  iterateCone(fInput, seedEta, seedPhi, 0, reduceConeSize, fOutput);
148 
149  }
150 
151 }
int i
Definition: DBlmapReader.cc:9
void iterateCone(const JetReco::InputCollection &fInput, double startRapidity, double startPhi, double startPt, bool reduceConeSize, InternalCollection *fOutput)
tuple cout
Definition: gather_cfg.py:121
void CMSMidpointAlgorithm::iterateCone ( const JetReco::InputCollection fInput,
double  startRapidity,
double  startPhi,
double  startPt,
bool  reduceConeSize,
InternalCollection fOutput 
)
private

Iterate the proto-jet center until it is stable. Called by findStableConesFromSeeds and findStableConesFromMidPoints but public method to allow testing and studies.

Definition at line 155 of file CMSMidpointAlgorithm.cc.

References funct::abs(), gather_cfg::cout, ProtoJet::e(), reco::Candidate::et(), reco::Candidate::eta(), ProtoJet::getTowerList(), i, ProtoJet::p4(), ProtoJet::phi(), reco::Candidate::phi(), ProtoJet::pt(), ProtoJet::putTowers(), mathSSE::sqrt(), lumiQTWidget::t, and ProtoJet::y().

160  {
161  // The workhorse of the algorithm.
162  // Throws a cone around a position and iterates the cone until the centroid and energy of the clsuter is stable.
163  // Uses a reduced cone size to prevent excessive cluster drift, ala CDF.
164  // Adds unique clusters to the list of proto-jets.
165  //
166  int nIterations = 0;
167  bool keepJet = true;
168  ProtoJet* trialCone = new ProtoJet;
169  double iterationtheConeRadius = theConeRadius;
170  if(reduceConeSize)iterationtheConeRadius *= sqrt(theConeAreaFraction);
171  while(++nIterations <= theMaxIterations + 1 && keepJet){
172 
173  // Last iteration uses the full cone size. Others use reduced cone size.
174  if(nIterations == theMaxIterations + 1)iterationtheConeRadius = theConeRadius;
175 
176  //Add all towers in cone and over threshold to the cluster
177  vector<InputItem> towersInSeedCluster
178  = towersWithinCone(fInput, startRapidity, startPhi, iterationtheConeRadius);
179  if(theDebugLevel>=2)cout << "[CMSMidpointAlgorithm] iter=" << nIterations << ", towers=" <<towersInSeedCluster.size();
180 
181  if(towersInSeedCluster.size()<1) { // Just in case there is an empty cone.
182  keepJet = false;
183  if(theDebugLevel>=2)cout << endl;
184  } else{
185  //Put seed cluster into trial cone
186  trialCone->putTowers(towersInSeedCluster);
187  double endRapidity = trialCone->y();
188  double endPhi = trialCone->phi();
189  double endPt = trialCone->pt();
190  double endE = trialCone->e();
191  if(theDebugLevel>=2)cout << ", y=" << endRapidity << ", phi=" << endPhi << ", PT=" << endPt << ", constituents=" << trialCone->getTowerList().size () << endl;
192  if(nIterations <= theMaxIterations){
193  // Do we have a stable cone?
194  if(std::abs(endRapidity-startRapidity)<.001 &&
195  std::abs(endPhi-startPhi)<.001 &&
196  std::abs(endE - startE) < .001) {
197  nIterations = theMaxIterations; // If cone size is reduced, then still one more iteration.
198  if(!reduceConeSize) ++nIterations; // Otherwise, this is the last iteration.
199  }
200  else{
201  // Another iteration.
202  startRapidity = endRapidity;
203  startPhi = endPhi;
204  startE = endE;
205  }
206  }
207  if(nIterations==theMaxIterations+1) {
208  for(vector<InputItem>::const_iterator i = towersInSeedCluster.begin(); i != towersInSeedCluster.end(); ++i) {
209  InputItem t = *i;
210  if(theDebugLevel>=2) cout << "[CMSMidpointAlgorithm] Tower " <<
211  i-towersInSeedCluster.begin() << ": eta=" << t->eta() <<
212  ", phi=" << t->phi() << ", ET=" << t->et() << endl;
213  }
214  }
215  }
216  }
217 
218  if(keepJet){ // Our trial cone is now stable
219  bool identical = false;
220  // Loop over proto-jets and check that our trial cone is a unique proto-jet
221  for (unsigned icone = 0; icone < stableCones->size(); icone++) {
222  if (trialCone->p4() == (*stableCones) [icone]->p4()) identical = true; // This proto-jet is not unique.
223  }
224  if(!identical){
225  stableCones->push_back(trialCone); // Save the unique proto-jets
226  trialCone = 0;
227  if(theDebugLevel>=2)cout << "[CMSMidpointAlgorithm] Unique Proto-Jet Saved" << endl;
228  }
229  }
230  delete trialCone;
231 }
int i
Definition: DBlmapReader.cc:9
double e() const
Definition: ProtoJet.h:67
virtual double et() const =0
transverse energy
virtual float eta() const =0
momentum pseudorapidity
const Constituents & getTowerList()
Definition: ProtoJet.cc:82
virtual float phi() const =0
momentum azimuthal angle
Transient Jet class used by the reconstruction algorithms.
Definition: ProtoJet.h:25
double y() const
Definition: ProtoJet.h:85
double pt() const
Definition: ProtoJet.h:75
double phi() const
Definition: ProtoJet.h:81
T sqrt(T t)
Definition: SSEVec.h:48
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const LorentzVector & p4() const
Definition: ProtoJet.h:91
tuple cout
Definition: gather_cfg.py:121
void putTowers(const Constituents &towers)
Definition: ProtoJet.cc:99
void CMSMidpointAlgorithm::run ( const JetReco::InputCollection fInput,
JetReco::OutputCollection fOutput 
)

Runs the algorithm and returns a list of caloJets. The user declares the vector and calls this method.

Definition at line 92 of file CMSMidpointAlgorithm.cc.

References dtNoiseDBValidation_cfg::cerr, gather_cfg::cout, and getHLTprescales::index.

93 {
94  if (!fOutput) {
95  std::cerr << "CMSMidpointAlgorithm::run-> ERROR: no output collection" << std::endl;
96  }
97  if(theDebugLevel>=1)cout << "[CMSMidpointAlgorithm] Num of input constituents = " << fInput.size() << endl;
98  if(theDebugLevel>=2) {
99  unsigned index = 0;
100  for (; index < fInput.size (); index++) {
101  cout << index << " consituent p/eta/phi: "
102  << fInput[index]->p() << '/'
103  << fInput[index]->eta() << '/'
104  << fInput[index]->phi() << std::endl;
105  }
106  }
107  // Find proto-jets from the seeds.
108  InternalCollection protoJets; // Initialize working container
109  // vector<ProtoJet> finalJets; // Final proto-jets container
110  findStableConesFromSeeds(fInput, &protoJets); //Find the proto-jets from the seeds
111  if(theDebugLevel>=1)cout << "[CMSMidpointAlgorithm] Num proto-jets from Seeds = " << protoJets.size() << endl;
112 
113  // Find proto-jets from the midpoints, and add them to the list
114  if(protoJets.size()>0)findStableConesFromMidPoints(fInput, &protoJets);// Add midpoints
115  if(theDebugLevel>=1)cout << "[CMSMidpointAlgorithm] Num proto-jets from Seeds and Midpoints = " << protoJets.size() << endl;
116 
117  // Split and merge the proto-jets, assigning each tower in the protojets to one and only one final jet.
118  if(protoJets.size()>0)splitAndMerge(fInput, &protoJets, fOutput); // Split and merge
119  if(theDebugLevel>=1)cout << "[CMSMidpointAlgorithm] Num final jets = " << fOutput->size() << endl;
120 
121  // Make the CaloJets from the final protojets
122  // MakeCaloJet(*theCtcp, finalJets, caloJets);
123 
124 }
void splitAndMerge(const JetReco::InputCollection &fInput, InternalCollection *fProtoJets, JetReco::OutputCollection *fFinalJets)
tuple cout
Definition: gather_cfg.py:121
void findStableConesFromSeeds(const JetReco::InputCollection &fInput, InternalCollection *fOutput)
void findStableConesFromMidPoints(const JetReco::InputCollection &fInput, InternalCollection *fOutput)
std::vector< ProtoJet * > InternalCollection
void CMSMidpointAlgorithm::splitAndMerge ( const JetReco::InputCollection fInput,
InternalCollection fProtoJets,
JetReco::OutputCollection fFinalJets 
)
private

Split and merge the proto-jets, assigning each tower in the protojets to one and only one final jet. Called by run but public method to allow testing and studies.

Definition at line 311 of file CMSMidpointAlgorithm.cc.

References dtNoiseDBValidation_cfg::cerr, gather_cfg::cout, Geom::deltaR2(), spr::find(), ProtoJet::getTowerList(), i, j, muon::overlap(), ProtoJet::p4(), ProtoJet::phi(), ProtoJet::pt(), ProtoJet::putTowers(), python.multivaluedict::sort(), and ProtoJet::y().

312  {
313  //
314  // This can use quite a bit of optimization and simplification.
315  //
316 
317  // Debugging
318  if(theDebugLevel>=2){
319  // Take a look at the proto-jets we were given
320  int numProtojets = stableCones->size();
321  for(int i = 0; i < numProtojets ; ++i){
322  const ProtoJet* icone = (*stableCones)[i];
323  int numTowers = icone->getTowerList().size();
324  cout << endl << "[CMSMidpointAlgorithm] ProtoJet " << i << ": PT=" << (*stableCones)[i]->pt()
325  << ", y="<< icone->y()
326  << ", phi="<< icone->phi()
327  << ", ntow="<< numTowers << endl;
328  ProtoJet::Constituents protojetTowers = icone->getTowerList();
329  for(int j = 0; j < numTowers; ++j){
330  cout << "[CMSMidpointAlgorithm] Tower " << j << ": ET=" << protojetTowers[j]->et()
331  << ", eta="<< protojetTowers[j]->eta()
332  << ", phi="<< protojetTowers[j]->phi() << endl;
333  }
334  }
335  }
336 
337  // Start of split and merge algorithm
338  bool mergingNotFinished = true;
339  while(mergingNotFinished){
340 
341  // Sort the stable cones (highest pt first).
342  GreaterByPtPtr<ProtoJet> compJets;
343  sort(stableCones->begin(),stableCones->end(),compJets);
344  // clean removed clusters
345  //InternalCollection::const_iterator i = find (stableCones->begin(), stableCones->end(), (ProtoJet*)0);
346  // for (; i != stableCones->end(); ++i) std::cout << "CMSMidpointAlgorithm::splitAndMerge-> removing pointer " << *i << std::endl;
347  stableCones->erase (find (stableCones->begin(), stableCones->end(), (ProtoJet*)0), stableCones->end());
348 
349  // Start with the highest pt cone left in list. List changes in loop,
350  // getting smaller with each iteration.
351  InternalCollection::iterator stableConeIter1 = stableCones->begin();
352 
353  if(stableConeIter1 == stableCones->end()) { // Stable cone list empty?
354  mergingNotFinished = false;
355  }
356  else {
357  ProtoJet* stableCone1 = *stableConeIter1;
358  if (!stableCone1) {
359  std::cerr << "CMSMidpointAlgorithm::splitAndMerge-> Error: stableCone1 should never be 0" << std::endl;
360  continue;
361  }
362  bool coneNotModified = true;
363 
364  // Determine whether highest pt cone has an overlap with other stable cones.
365  InternalCollection::iterator stableConeIter2 = stableConeIter1;
366  ++stableConeIter2; // Iterator for 2nd highest pt cone, just past 1st cone.
367 
368  while(coneNotModified && stableConeIter2 != stableCones->end()){
369  ProtoJet* stableCone2 = *stableConeIter2;
370  if (!stableCone2) {
371  std::cerr << "CMSMidpointAlgorithm::splitAndMerge-> Error: stableCone2 should never be 0" << std::endl;
372  continue;
373  }
374 
375  // Calculate overlap of the two cones.
376  bool overlap = false;
377  vector<InputItem> overlapTowers; // Make a list to hold the overlap towers
378  //cout << "1st cone num towers=" << stableCone1->getTowerList().size() << endl;
379  //int numTowers1=0;
380 
381  //Loop over towers in higher Pt cone
382  for(ProtoJet::Constituents::const_iterator towerIter1 = stableCone1->getTowerList().begin();
383  towerIter1 != stableCone1->getTowerList().end();
384  ++towerIter1){
385  //cout << "1st cone tower " << numTowers1 << endl;
386  //++numTowers1;
387  bool isInCone2 = false;
388  //cout << "2nd cone num towers=" << stableCone2->getTowerList().size() << endl;
389  //int numTowers2=0;
390 
391  // Loop over towers in lower Pt cone
392  for(ProtoJet::Constituents::const_iterator towerIter2 = stableCone2->getTowerList().begin();
393  towerIter2 != stableCone2->getTowerList().end();
394  ++towerIter2) {
395  //cout << "2st cone tower " << numTowers2 << endl;
396  //++numTowers2;
397 
398  // Check if towers are the same by checking for unique eta, phi and energy values.
399  // Will want to replace this with checking the tower index when available.
400  if(sameTower(*towerIter1, *towerIter2)) {
401  isInCone2 = true; //The tower is in both cones
402  //cout << "Merging found overlap tower: eta=" << (*towerIter1)->eta << ", phi=" << (*towerIter1)->phi << " in 1st protojet " <<
403  // " and eta=" << (*towerIter2)->eta << ", phi=" << (*towerIter2)->phi << " in 2nd protojet " << endl;
404  }
405  }
406  if(isInCone2){
407  overlapTowers.push_back(*towerIter1); //Add tower in both cones to the overlap list
408  overlap = true;
409  }
410  }
411  if(overlap){
412  // non-empty overlap. Decide on splitting or merging.
413 
414  // Make a proto-jet with the overlap towers so we can calculate things for the overlap
416  overlap.putTowers(overlapTowers);
417  coneNotModified = false;
418 
419  // Compare the overlap pt with the overlap fractcion threshold times the lower jet pt.
420  if(overlap.pt() >= theOverlapThreshold*stableCone2->pt()){
421 
422  // Merge the two cones.
423  // Get a copy of the list of towers in higher Pt proto-jet
424  ProtoJet::Constituents stableCone1Towers = stableCone1->getTowerList();
425 
426  //Loop over the list of towers lower Pt jet
427  for(ProtoJet::Constituents::const_iterator towerIter2 = stableCone2->getTowerList().begin();
428  towerIter2 != stableCone2->getTowerList().end();
429  ++towerIter2){
430  bool isInOverlap = false;
431 
432  //Check if that tower is in the overlap region
433  for(vector<InputItem>::iterator overlapTowerIter = overlapTowers.begin();
434  overlapTowerIter != overlapTowers.end();
435  ++overlapTowerIter){
436  // Check if towers are the same by checking for unique eta, phi and energy values.
437  // Will want to replace this with checking the tower index when available.
438  if(sameTower(*overlapTowerIter, *towerIter2)) isInOverlap = true;
439  }
440  //Add non-overlap tower from proto-jet 2 into the proto-jet 1 tower list.
441  if(!isInOverlap)stableCone1Towers.push_back(*towerIter2);
442  }
443 
444  if(theDebugLevel>=2)cout << endl << "[CMSMidpointAlgorithm] Merging: 1st Proto-jet grows: "
445  " y=" << stableCone1->y() <<
446  ", phi=" << stableCone1->phi() <<
447  " increases from " << stableCone1->getTowerList().size() <<
448  " to " << stableCone1Towers.size() << " towers." << endl;
449 
450  // Put the new expanded list of towers into the first proto-jet
451  stableCone1->putTowers(stableCone1Towers);
452 
453  if(theDebugLevel>=2)cout << "[CMSMidpointAlgorithm] Merging: 1st protojet now at y=" << stableCone1->y() <<
454  ", phi=" << stableCone1->phi() << endl;
455 
456  if(theDebugLevel>=2)cout << "[CMSMidpointAlgorithm] Merging: 2nd Proto-jet removed:"
457  " y=" << stableCone2->y() <<
458  ", phi=" << stableCone2->phi() << endl;
459 
460  // Remove the second proto-jet.
461  delete *stableConeIter2;
462  *stableConeIter2 = 0;
463 
464  }
465  else{
466  // Split the two proto-jets.
467 
468  // Create lists of towers to remove from each proto-jet
469  vector<InputItem> removeFromCone1,removeFromCone2;
470 
471  // Which tower goes where?
472  // Loop over the overlap towers
473  for(vector<InputItem>::iterator towerIter = overlapTowers.begin();
474  towerIter != overlapTowers.end();
475  ++towerIter){
476  double dR2Jet1 = deltaR2 ((*towerIter)->p4().Rapidity(), (*towerIter)->phi(),
477  stableCone1->y(), stableCone1->phi());
478  // Calculate distance from proto-jet 2.
479  double dR2Jet2 = deltaR2 ((*towerIter)->p4().Rapidity(), (*towerIter)->phi(),
480  stableCone2->y(), stableCone2->phi());
481 
482  if(dR2Jet1 < dR2Jet2){
483  // Tower is closer to proto-jet 1. To be removed from proto-jet 2.
484  removeFromCone2.push_back(*towerIter);
485  }
486  else {
487  // Tower is closer to proto-jet 2. To be removed from proto-jet 1.
488  removeFromCone1.push_back(*towerIter);
489  }
490  }
491  // Remove towers in the overlap region from the cones to which they have the larger distance.
492 
493  // Remove towers from proto-jet 1.
494  vector<InputItem> towerList1 (stableCone1->getTowerList().begin(), stableCone1->getTowerList().end());
495 
496  // Loop over towers in remove list
497  for(vector<InputItem>::iterator towerIter = removeFromCone1.begin();
498  towerIter != removeFromCone1.end();
499  ++towerIter) {
500  // Loop over towers in protojet
501  for(vector<InputItem>::iterator towerIter1 = towerList1.begin(); towerIter1 != towerList1.end(); ++towerIter1) {
502  // Check if they are equal
503  if(sameTower(*towerIter, *towerIter1)) {
504  // Remove the tower
505  towerList1.erase(towerIter1);
506  break;
507  }
508  }
509  }
510 
511  if(theDebugLevel>=2)cout << endl << "[CMSMidpointAlgorithm] Splitting: 1st Proto-jet shrinks: y=" <<
512  stableCone1->y() <<
513  ", phi=" << stableCone1->phi() <<
514  " decreases from" << stableCone1->getTowerList().size() <<
515  " to " << towerList1.size() << " towers." << endl;
516 
517  //Put the new reduced list of towers into proto-jet 1.
518  stableCone1->putTowers(towerList1);
519  // Remove towers from cone 2.
520  vector<InputItem> towerList2 (stableCone2->getTowerList().begin(), stableCone2->getTowerList().end());
521 
522  // Loop over towers in remove list
523  for(vector<InputItem>::iterator towerIter = removeFromCone2.begin();
524  towerIter != removeFromCone2.end();
525  ++towerIter) {
526  // Loop over towers in protojet
527  for(vector<InputItem>::iterator towerIter2 = towerList2.begin(); towerIter2 != towerList2.end(); ++towerIter2){
528  // Check if they are equal
529  if(sameTower(*towerIter, *towerIter2)) {
530  // Remove the tower
531  towerList2.erase(towerIter2);
532  break;
533  }
534  }
535  }
536 
537  if(theDebugLevel>=2)cout << "[CMSMidpointAlgorithm] Splitting: 2nd Proto-jet shrinks: y=" <<
538  stableCone2->y() <<
539  ", phi=" << stableCone2->phi() <<
540  " decreases from" << stableCone2->getTowerList().size() <<
541  " to " << towerList2.size() << " towers." << endl;
542 
543  //Put the new reduced list of towers into proto-jet 2.
544  stableCone2->putTowers(towerList2);
545  }
546  }
547  else {
548  if(theDebugLevel>=2)cout << endl <<
549  "[CMSMidpointAlgorithm] no overlap between 1st protojet at y=" << stableCone1->y() <<
550  ", phi=" << stableCone1->phi() <<
551  " and 2nd protojet at y=" << stableCone2->y() <<
552  ", phi=" << stableCone2->phi() <<endl;
553  }
554 
555  ++stableConeIter2; //Increment iterator to the next highest Pt protojet
556  }
557  if(coneNotModified){
558 
559  if(theDebugLevel>=2)cout <<
560  "[CMSMidpointAlgorithm] Saving: Proto-jet at y=" << stableCone1->y() <<
561  ", phi=" << stableCone1->phi() << " has no overlap" << endl;
562 
563  fFinalJets->push_back(ProtoJet (stableCone1->p4(), stableCone1->getTowerList()));
564  delete *stableConeIter1;
565  *stableConeIter1 = 0;
566  }
567  }
568  }
569 
570  GreaterByPt<ProtoJet> compJets;
571  sort(fFinalJets->begin(),fFinalJets->end(),compJets);
572 }
int i
Definition: DBlmapReader.cc:9
JetReco::InputCollection Constituents
Definition: ProtoJet.h:28
const Constituents & getTowerList()
Definition: ProtoJet.cc:82
Transient Jet class used by the reconstruction algorithms.
Definition: ProtoJet.h:25
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
double y() const
Definition: ProtoJet.h:85
bool overlap(const reco::Muon &muon1, const reco::Muon &muon2, double pullX=1.0, double pullY=1.0, bool checkAdjacentChambers=false)
double pt() const
Definition: ProtoJet.h:75
double phi() const
Definition: ProtoJet.h:81
int j
Definition: DBlmapReader.cc:9
double deltaR2(const Vector1 &v1, const Vector2 &v2)
Definition: VectorUtil.h:78
const LorentzVector & p4() const
Definition: ProtoJet.h:91
tuple cout
Definition: gather_cfg.py:121
void putTowers(const Constituents &towers)
Definition: ProtoJet.cc:99

Member Data Documentation

double CMSMidpointAlgorithm::theConeAreaFraction
private

Definition at line 108 of file CMSMidpointAlgorithm.h.

double CMSMidpointAlgorithm::theConeRadius
private

Definition at line 107 of file CMSMidpointAlgorithm.h.

int CMSMidpointAlgorithm::theDebugLevel
private

Definition at line 112 of file CMSMidpointAlgorithm.h.

int CMSMidpointAlgorithm::theMaxIterations
private

Definition at line 110 of file CMSMidpointAlgorithm.h.

int CMSMidpointAlgorithm::theMaxPairSize
private

Definition at line 109 of file CMSMidpointAlgorithm.h.

double CMSMidpointAlgorithm::theOverlapThreshold
private

Definition at line 111 of file CMSMidpointAlgorithm.h.

double CMSMidpointAlgorithm::theSeedThreshold
private

Definition at line 106 of file CMSMidpointAlgorithm.h.