00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034 #include "DataFormats/Candidate/interface/Candidate.h"
00035 #include "DataFormats/Math/interface/deltaR.h"
00036 #include "RecoParticleFlow/PFRootEvent/interface/ProtoJet.h"
00037 #include "RecoParticleFlow/PFRootEvent/interface/JetAlgoHelper.h"
00038 #include "RecoParticleFlow/PFRootEvent/interface/CMSMidpointAlgorithm.h"
00039
00040
00041
00042 using namespace std;
00043 using namespace reco;
00044 using namespace JetReco;
00045
00046
00047 namespace {
00048
00049 bool sameTower (InputItem c1, InputItem c2) {
00050 return c1 == c2;
00051 }
00052
00053 std::vector<InputItem> towersWithinCone(const InputCollection& fInput, double coneEta, double conePhi, double coneRadius){
00054 std::vector<InputItem> result;
00055 double r2 = coneRadius*coneRadius;
00056 InputCollection::const_iterator towerIter = fInput.begin();
00057 InputCollection::const_iterator towerIterEnd = fInput.end();
00058 for (;towerIter != towerIterEnd; ++towerIter) {
00059 InputItem caloTowerPointer = *towerIter;
00060 double dR2 = deltaR2 (coneEta, conePhi, caloTowerPointer->eta(), caloTowerPointer->phi());
00061 if(dR2 < r2){
00062 result.push_back(caloTowerPointer);
00063 }
00064 }
00065 return result;
00066 }
00067
00068
00069 std::vector<InputItem> etOrderedCaloTowers(const InputCollection& fInput, double etThreshold) {
00070 std::vector<InputItem> result;
00071 InputCollection::const_iterator towerIter = fInput.begin();
00072 InputCollection::const_iterator towerIterEnd = fInput.end();
00073 for (;towerIter != towerIterEnd; ++towerIter) {
00074 InputItem caloTowerPointer = *towerIter;
00075 if(caloTowerPointer->et() > etThreshold){
00076 result.push_back(caloTowerPointer);
00077 }
00078 }
00079 GreaterByEtRef <InputItem> compCandidate;
00080 sort (result.begin(), result.end(), compCandidate);
00081 return result;
00082 }
00083
00084
00085
00086
00087 }
00088
00089
00090
00091
00092 void CMSMidpointAlgorithm::run(const InputCollection& fInput, OutputCollection* fOutput)
00093 {
00094 if (!fOutput) {
00095 std::cerr << "CMSMidpointAlgorithm::run-> ERROR: no output collection" << std::endl;
00096 }
00097 if(theDebugLevel>=1)cout << "[CMSMidpointAlgorithm] Num of input constituents = " << fInput.size() << endl;
00098 if(theDebugLevel>=2) {
00099 unsigned index = 0;
00100 for (; index < fInput.size (); index++) {
00101 cout << index << " consituent p/eta/phi: "
00102 << fInput[index]->p() << '/'
00103 << fInput[index]->eta() << '/'
00104 << fInput[index]->phi() << std::endl;
00105 }
00106 }
00107
00108 InternalCollection protoJets;
00109
00110 findStableConesFromSeeds(fInput, &protoJets);
00111 if(theDebugLevel>=1)cout << "[CMSMidpointAlgorithm] Num proto-jets from Seeds = " << protoJets.size() << endl;
00112
00113
00114 if(protoJets.size()>0)findStableConesFromMidPoints(fInput, &protoJets);
00115 if(theDebugLevel>=1)cout << "[CMSMidpointAlgorithm] Num proto-jets from Seeds and Midpoints = " << protoJets.size() << endl;
00116
00117
00118 if(protoJets.size()>0)splitAndMerge(fInput, &protoJets, fOutput);
00119 if(theDebugLevel>=1)cout << "[CMSMidpointAlgorithm] Num final jets = " << fOutput->size() << endl;
00120
00121
00122
00123
00124 }
00125
00126
00127
00128
00129
00130 void CMSMidpointAlgorithm::findStableConesFromSeeds(const InputCollection& fInput, InternalCollection* fOutput) {
00131
00132
00133
00134 bool reduceConeSize = true;
00135
00136
00137 vector<InputItem> seedTowers = etOrderedCaloTowers(fInput, theSeedThreshold);
00138
00139
00140 for(vector<InputItem>::const_iterator i = seedTowers.begin(); i != seedTowers.end(); ++i) {
00141 double seedEta = (*i)->eta ();
00142 double seedPhi = (*i)->phi ();
00143
00144
00145
00146 if(theDebugLevel>=2) cout << endl << "[CMSMidpointAlgorithm] seed " << i-seedTowers.begin() << ": eta=" << seedEta << ", phi=" << seedPhi << endl;
00147 iterateCone(fInput, seedEta, seedPhi, 0, reduceConeSize, fOutput);
00148
00149 }
00150
00151 }
00152
00153
00154
00155 void CMSMidpointAlgorithm::iterateCone(const InputCollection& fInput,
00156 double startRapidity,
00157 double startPhi,
00158 double startE,
00159 bool reduceConeSize,
00160 InternalCollection* stableCones) {
00161
00162
00163
00164
00165
00166 int nIterations = 0;
00167 bool keepJet = true;
00168 ProtoJet* trialCone = new ProtoJet;
00169 double iterationtheConeRadius = theConeRadius;
00170 if(reduceConeSize)iterationtheConeRadius *= sqrt(theConeAreaFraction);
00171 while(++nIterations <= theMaxIterations + 1 && keepJet){
00172
00173
00174 if(nIterations == theMaxIterations + 1)iterationtheConeRadius = theConeRadius;
00175
00176
00177 vector<InputItem> towersInSeedCluster
00178 = towersWithinCone(fInput, startRapidity, startPhi, iterationtheConeRadius);
00179 if(theDebugLevel>=2)cout << "[CMSMidpointAlgorithm] iter=" << nIterations << ", towers=" <<towersInSeedCluster.size();
00180
00181 if(towersInSeedCluster.size()<1) {
00182 keepJet = false;
00183 if(theDebugLevel>=2)cout << endl;
00184 } else{
00185
00186 trialCone->putTowers(towersInSeedCluster);
00187 double endRapidity = trialCone->y();
00188 double endPhi = trialCone->phi();
00189 double endPt = trialCone->pt();
00190 double endE = trialCone->e();
00191 if(theDebugLevel>=2)cout << ", y=" << endRapidity << ", phi=" << endPhi << ", PT=" << endPt << ", constituents=" << trialCone->getTowerList().size () << endl;
00192 if(nIterations <= theMaxIterations){
00193
00194 if(std::abs(endRapidity-startRapidity)<.001 &&
00195 std::abs(endPhi-startPhi)<.001 &&
00196 std::abs(endE - startE) < .001) {
00197 nIterations = theMaxIterations;
00198 if(!reduceConeSize) ++nIterations;
00199 }
00200 else{
00201
00202 startRapidity = endRapidity;
00203 startPhi = endPhi;
00204 startE = endE;
00205 }
00206 }
00207 if(nIterations==theMaxIterations+1) {
00208 for(vector<InputItem>::const_iterator i = towersInSeedCluster.begin(); i != towersInSeedCluster.end(); ++i) {
00209 InputItem t = *i;
00210 if(theDebugLevel>=2) cout << "[CMSMidpointAlgorithm] Tower " <<
00211 i-towersInSeedCluster.begin() << ": eta=" << t->eta() <<
00212 ", phi=" << t->phi() << ", ET=" << t->et() << endl;
00213 }
00214 }
00215 }
00216 }
00217
00218 if(keepJet){
00219 bool identical = false;
00220
00221 for (unsigned icone = 0; icone < stableCones->size(); icone++) {
00222 if (trialCone->p4() == (*stableCones) [icone]->p4()) identical = true;
00223 }
00224 if(!identical){
00225 stableCones->push_back(trialCone);
00226 trialCone = 0;
00227 if(theDebugLevel>=2)cout << "[CMSMidpointAlgorithm] Unique Proto-Jet Saved" << endl;
00228 }
00229 }
00230 delete trialCone;
00231 }
00232
00233
00234
00235
00236 void CMSMidpointAlgorithm::findStableConesFromMidPoints(const InputCollection& fInput, InternalCollection* stableCones){
00237
00238
00239
00240 vector< vector<bool> > distanceOK;
00241 distanceOK.resize(stableCones->size() - 1);
00242 for(unsigned int nCluster1 = 1; nCluster1 < stableCones->size(); ++nCluster1){
00243 distanceOK[nCluster1 - 1].resize(nCluster1);
00244 const ProtoJet* cluster1 = (*stableCones)[nCluster1];
00245 for(unsigned int nCluster2 = 0; nCluster2 < nCluster1; ++nCluster2){
00246 const ProtoJet* cluster2 = (*stableCones)[nCluster2];
00247 double dR2 = deltaR2 (cluster1->y(), cluster1->phi(), cluster2->y(), cluster2->phi());
00248 distanceOK[nCluster1 - 1][nCluster2] = dR2 < 4*theConeRadius*theConeRadius;
00249 }
00250 }
00251
00252
00253 vector< vector<int> > pairs(0);
00254 vector<int> testPair(0);
00255 int maxClustersInPair = theMaxPairSize;
00256 if(!maxClustersInPair)maxClustersInPair = stableCones->size();
00257 addClustersToPairs(fInput, testPair,pairs,distanceOK,maxClustersInPair);
00258
00259
00260 bool reduceConeSize = false;
00261 for(vector<vector<int> >::const_iterator iPair = pairs.begin(); iPair != pairs.end(); ++iPair) {
00262 const vector<int> & Pair = *iPair;
00263
00264 reco::Particle::LorentzVector midPoint (0,0,0,0);
00265 for(vector<int>::const_iterator iPairMember = Pair.begin(); iPairMember != Pair.end(); ++iPairMember) {
00266 midPoint += (*stableCones)[*iPairMember]->p4();
00267 }
00268 if(theDebugLevel>=2)cout << endl << "[CMSMidpointAlgorithm] midpoint " << iPair-pairs.begin() << ": y = " << midPoint.Rapidity() << ", phi=" << midPoint.Phi() <<
00269 ", size=" << Pair.size() << endl;
00270 iterateCone(fInput, midPoint.Rapidity(),midPoint.Phi(),midPoint.e(),reduceConeSize,stableCones);
00271 }
00272 GreaterByPtPtr<ProtoJet> compJets;
00273 sort (stableCones->begin(), stableCones->end(), compJets);
00274 }
00275
00276
00277
00278 void CMSMidpointAlgorithm::addClustersToPairs(const InputCollection& fInput,
00279 vector<int>& testPair, vector< vector<int> >& pairs,
00280 vector< vector<bool> >& distanceOK, int maxClustersInPair)
00281 {
00282
00283
00284
00285 int nextClusterStart = 0;
00286 if(testPair.size())
00287 nextClusterStart = testPair.back() + 1;
00288 for(unsigned int nextCluster = nextClusterStart; nextCluster <= distanceOK.size(); ++nextCluster){
00289
00290 bool addCluster = true;
00291 for(unsigned int iCluster = 0; iCluster < testPair.size() && addCluster; ++iCluster)
00292 if(!distanceOK[nextCluster - 1][testPair[iCluster]])
00293 addCluster = false;
00294 if(addCluster){
00295
00296 testPair.push_back(nextCluster);
00297
00298 if(testPair.size() > 1)
00299 pairs.push_back(testPair);
00300
00301 if(testPair.size() < unsigned(maxClustersInPair))
00302 addClustersToPairs(fInput, testPair,pairs,distanceOK,maxClustersInPair);
00303
00304 testPair.pop_back();
00305 }
00306 }
00307 }
00308
00309
00310
00311 void CMSMidpointAlgorithm::splitAndMerge(const InputCollection& fInput,
00312 InternalCollection* stableCones, OutputCollection* fFinalJets) {
00313
00314
00315
00316
00317
00318 if(theDebugLevel>=2){
00319
00320 int numProtojets = stableCones->size();
00321 for(int i = 0; i < numProtojets ; ++i){
00322 const ProtoJet* icone = (*stableCones)[i];
00323 int numTowers = icone->getTowerList().size();
00324 cout << endl << "[CMSMidpointAlgorithm] ProtoJet " << i << ": PT=" << (*stableCones)[i]->pt()
00325 << ", y="<< icone->y()
00326 << ", phi="<< icone->phi()
00327 << ", ntow="<< numTowers << endl;
00328 ProtoJet::Constituents protojetTowers = icone->getTowerList();
00329 for(int j = 0; j < numTowers; ++j){
00330 cout << "[CMSMidpointAlgorithm] Tower " << j << ": ET=" << protojetTowers[j]->et()
00331 << ", eta="<< protojetTowers[j]->eta()
00332 << ", phi="<< protojetTowers[j]->phi() << endl;
00333 }
00334 }
00335 }
00336
00337
00338 bool mergingNotFinished = true;
00339 while(mergingNotFinished){
00340
00341
00342 GreaterByPtPtr<ProtoJet> compJets;
00343 sort(stableCones->begin(),stableCones->end(),compJets);
00344
00345
00346
00347 stableCones->erase (find (stableCones->begin(), stableCones->end(), (ProtoJet*)0), stableCones->end());
00348
00349
00350
00351 InternalCollection::iterator stableConeIter1 = stableCones->begin();
00352
00353 if(stableConeIter1 == stableCones->end()) {
00354 mergingNotFinished = false;
00355 }
00356 else {
00357 ProtoJet* stableCone1 = *stableConeIter1;
00358 if (!stableCone1) {
00359 std::cerr << "CMSMidpointAlgorithm::splitAndMerge-> Error: stableCone1 should never be 0" << std::endl;
00360 continue;
00361 }
00362 bool coneNotModified = true;
00363
00364
00365 InternalCollection::iterator stableConeIter2 = stableConeIter1;
00366 ++stableConeIter2;
00367
00368 while(coneNotModified && stableConeIter2 != stableCones->end()){
00369 ProtoJet* stableCone2 = *stableConeIter2;
00370 if (!stableCone2) {
00371 std::cerr << "CMSMidpointAlgorithm::splitAndMerge-> Error: stableCone2 should never be 0" << std::endl;
00372 continue;
00373 }
00374
00375
00376 bool overlap = false;
00377 vector<InputItem> overlapTowers;
00378
00379
00380
00381
00382 for(ProtoJet::Constituents::const_iterator towerIter1 = stableCone1->getTowerList().begin();
00383 towerIter1 != stableCone1->getTowerList().end();
00384 ++towerIter1){
00385
00386
00387 bool isInCone2 = false;
00388
00389
00390
00391
00392 for(ProtoJet::Constituents::const_iterator towerIter2 = stableCone2->getTowerList().begin();
00393 towerIter2 != stableCone2->getTowerList().end();
00394 ++towerIter2) {
00395
00396
00397
00398
00399
00400 if(sameTower(*towerIter1, *towerIter2)) {
00401 isInCone2 = true;
00402
00403
00404 }
00405 }
00406 if(isInCone2){
00407 overlapTowers.push_back(*towerIter1);
00408 overlap = true;
00409 }
00410 }
00411 if(overlap){
00412
00413
00414
00415 ProtoJet overlap;
00416 overlap.putTowers(overlapTowers);
00417 coneNotModified = false;
00418
00419
00420 if(overlap.pt() >= theOverlapThreshold*stableCone2->pt()){
00421
00422
00423
00424 ProtoJet::Constituents stableCone1Towers = stableCone1->getTowerList();
00425
00426
00427 for(ProtoJet::Constituents::const_iterator towerIter2 = stableCone2->getTowerList().begin();
00428 towerIter2 != stableCone2->getTowerList().end();
00429 ++towerIter2){
00430 bool isInOverlap = false;
00431
00432
00433 for(vector<InputItem>::iterator overlapTowerIter = overlapTowers.begin();
00434 overlapTowerIter != overlapTowers.end();
00435 ++overlapTowerIter){
00436
00437
00438 if(sameTower(*overlapTowerIter, *towerIter2)) isInOverlap = true;
00439 }
00440
00441 if(!isInOverlap)stableCone1Towers.push_back(*towerIter2);
00442 }
00443
00444 if(theDebugLevel>=2)cout << endl << "[CMSMidpointAlgorithm] Merging: 1st Proto-jet grows: "
00445 " y=" << stableCone1->y() <<
00446 ", phi=" << stableCone1->phi() <<
00447 " increases from " << stableCone1->getTowerList().size() <<
00448 " to " << stableCone1Towers.size() << " towers." << endl;
00449
00450
00451 stableCone1->putTowers(stableCone1Towers);
00452
00453 if(theDebugLevel>=2)cout << "[CMSMidpointAlgorithm] Merging: 1st protojet now at y=" << stableCone1->y() <<
00454 ", phi=" << stableCone1->phi() << endl;
00455
00456 if(theDebugLevel>=2)cout << "[CMSMidpointAlgorithm] Merging: 2nd Proto-jet removed:"
00457 " y=" << stableCone2->y() <<
00458 ", phi=" << stableCone2->phi() << endl;
00459
00460
00461 delete *stableConeIter2;
00462 *stableConeIter2 = 0;
00463
00464 }
00465 else{
00466
00467
00468
00469 vector<InputItem> removeFromCone1,removeFromCone2;
00470
00471
00472
00473 for(vector<InputItem>::iterator towerIter = overlapTowers.begin();
00474 towerIter != overlapTowers.end();
00475 ++towerIter){
00476 double dR2Jet1 = deltaR2 ((*towerIter)->p4().Rapidity(), (*towerIter)->phi(),
00477 stableCone1->y(), stableCone1->phi());
00478
00479 double dR2Jet2 = deltaR2 ((*towerIter)->p4().Rapidity(), (*towerIter)->phi(),
00480 stableCone2->y(), stableCone2->phi());
00481
00482 if(dR2Jet1 < dR2Jet2){
00483
00484 removeFromCone2.push_back(*towerIter);
00485 }
00486 else {
00487
00488 removeFromCone1.push_back(*towerIter);
00489 }
00490 }
00491
00492
00493
00494 vector<InputItem> towerList1 (stableCone1->getTowerList().begin(), stableCone1->getTowerList().end());
00495
00496
00497 for(vector<InputItem>::iterator towerIter = removeFromCone1.begin();
00498 towerIter != removeFromCone1.end();
00499 ++towerIter) {
00500
00501 for(vector<InputItem>::iterator towerIter1 = towerList1.begin(); towerIter1 != towerList1.end(); ++towerIter1) {
00502
00503 if(sameTower(*towerIter, *towerIter1)) {
00504
00505 towerList1.erase(towerIter1);
00506 break;
00507 }
00508 }
00509 }
00510
00511 if(theDebugLevel>=2)cout << endl << "[CMSMidpointAlgorithm] Splitting: 1st Proto-jet shrinks: y=" <<
00512 stableCone1->y() <<
00513 ", phi=" << stableCone1->phi() <<
00514 " decreases from" << stableCone1->getTowerList().size() <<
00515 " to " << towerList1.size() << " towers." << endl;
00516
00517
00518 stableCone1->putTowers(towerList1);
00519
00520 vector<InputItem> towerList2 (stableCone2->getTowerList().begin(), stableCone2->getTowerList().end());
00521
00522
00523 for(vector<InputItem>::iterator towerIter = removeFromCone2.begin();
00524 towerIter != removeFromCone2.end();
00525 ++towerIter) {
00526
00527 for(vector<InputItem>::iterator towerIter2 = towerList2.begin(); towerIter2 != towerList2.end(); ++towerIter2){
00528
00529 if(sameTower(*towerIter, *towerIter2)) {
00530
00531 towerList2.erase(towerIter2);
00532 break;
00533 }
00534 }
00535 }
00536
00537 if(theDebugLevel>=2)cout << "[CMSMidpointAlgorithm] Splitting: 2nd Proto-jet shrinks: y=" <<
00538 stableCone2->y() <<
00539 ", phi=" << stableCone2->phi() <<
00540 " decreases from" << stableCone2->getTowerList().size() <<
00541 " to " << towerList2.size() << " towers." << endl;
00542
00543
00544 stableCone2->putTowers(towerList2);
00545 }
00546 }
00547 else {
00548 if(theDebugLevel>=2)cout << endl <<
00549 "[CMSMidpointAlgorithm] no overlap between 1st protojet at y=" << stableCone1->y() <<
00550 ", phi=" << stableCone1->phi() <<
00551 " and 2nd protojet at y=" << stableCone2->y() <<
00552 ", phi=" << stableCone2->phi() <<endl;
00553 }
00554
00555 ++stableConeIter2;
00556 }
00557 if(coneNotModified){
00558
00559 if(theDebugLevel>=2)cout <<
00560 "[CMSMidpointAlgorithm] Saving: Proto-jet at y=" << stableCone1->y() <<
00561 ", phi=" << stableCone1->phi() << " has no overlap" << endl;
00562
00563 fFinalJets->push_back(ProtoJet (stableCone1->p4(), stableCone1->getTowerList()));
00564 delete *stableConeIter1;
00565 *stableConeIter1 = 0;
00566 }
00567 }
00568 }
00569
00570 GreaterByPt<ProtoJet> compJets;
00571 sort(fFinalJets->begin(),fFinalJets->end(),compJets);
00572 }