00001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00002
00003 #include "RecoParticleFlow/PFProducer/interface/PFAlgo.h"
00004 #include "RecoParticleFlow/PFProducer/interface/PFMuonAlgo.h"
00005 #include "RecoParticleFlow/PFProducer/interface/PFElectronAlgo.h"
00006 #include "RecoParticleFlow/PFProducer/interface/PFPhotonAlgo.h"
00007 #include "RecoParticleFlow/PFProducer/interface/PFElectronExtraEqual.h"
00008
00009 #include "RecoParticleFlow/PFClusterTools/interface/PFEnergyCalibration.h"
00010 #include "RecoParticleFlow/PFClusterTools/interface/PFEnergyCalibrationHF.h"
00011 #include "RecoParticleFlow/PFClusterTools/interface/PFSCEnergyCalibration.h"
00012
00013 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
00014 #include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
00015 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementTrack.h"
00016 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementCluster.h"
00017 #include "DataFormats/ParticleFlowReco/interface/PFRecTrack.h"
00018 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
00019 #include "DataFormats/ParticleFlowReco/interface/PFLayer.h"
00020
00021 #include "DataFormats/VertexReco/interface/Vertex.h"
00022 #include "DataFormats/TrackReco/interface/Track.h"
00023 #include "DataFormats/MuonReco/interface/Muon.h"
00024
00025
00026 #include "DataFormats/Common/interface/OrphanHandle.h"
00027
00028 #include "DataFormats/Provenance/interface/ProductID.h"
00029 #include "DataFormats/Math/interface/LorentzVector.h"
00030
00031
00032 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00033 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
00034
00035 #include "DataFormats/ParticleFlowReco/interface/PFDisplacedVertex.h"
00036 #include "DataFormats/ParticleFlowReco/interface/PFDisplacedVertexFwd.h"
00037
00038
00039 #include "Math/PxPyPzM4D.h"
00040 #include "Math/LorentzVector.h"
00041 #include "Math/DisplacementVector3D.h"
00042 #include "Math/SMatrix.h"
00043 #include "TDecompChol.h"
00044
00045 #include "boost/graph/adjacency_matrix.hpp"
00046 #include "boost/graph/graph_utility.hpp"
00047
00048
00049 using namespace std;
00050 using namespace reco;
00051 using namespace boost;
00052
00053 typedef std::list< reco::PFBlockRef >::iterator IBR;
00054
00055
00056
00057 PFAlgo::PFAlgo()
00058 : pfCandidates_( new PFCandidateCollection),
00059 nSigmaECAL_(0),
00060 nSigmaHCAL_(1),
00061 algo_(1),
00062 debug_(false),
00063 pfele_(0),
00064 pfpho_(0),
00065 useVertices_(false)
00066 {}
00067
00068 PFAlgo::~PFAlgo() {
00069 if (usePFElectrons_) delete pfele_;
00070 if (usePFPhotons_) delete pfpho_;
00071 }
00072
00073
00074 void
00075 PFAlgo::setParameters(double nSigmaECAL,
00076 double nSigmaHCAL,
00077 const boost::shared_ptr<PFEnergyCalibration>& calibration,
00078 const boost::shared_ptr<PFEnergyCalibrationHF>& thepfEnergyCalibrationHF) {
00079
00080 nSigmaECAL_ = nSigmaECAL;
00081 nSigmaHCAL_ = nSigmaHCAL;
00082
00083 calibration_ = calibration;
00084 thepfEnergyCalibrationHF_ = thepfEnergyCalibrationHF;
00085
00086 }
00087
00088
00089 PFMuonAlgo* PFAlgo::getPFMuonAlgo() {
00090 return pfmu_;
00091 }
00092
00093
00094 void
00095 PFAlgo::setPFEleParameters(double mvaEleCut,
00096 string mvaWeightFileEleID,
00097 bool usePFElectrons,
00098 const boost::shared_ptr<PFSCEnergyCalibration>& thePFSCEnergyCalibration,
00099 const boost::shared_ptr<PFEnergyCalibration>& thePFEnergyCalibration,
00100 double sumEtEcalIsoForEgammaSC_barrel,
00101 double sumEtEcalIsoForEgammaSC_endcap,
00102 double coneEcalIsoForEgammaSC,
00103 double sumPtTrackIsoForEgammaSC_barrel,
00104 double sumPtTrackIsoForEgammaSC_endcap,
00105 unsigned int nTrackIsoForEgammaSC,
00106 double coneTrackIsoForEgammaSC,
00107 bool applyCrackCorrections,
00108 bool usePFSCEleCalib,
00109 bool useEGElectrons,
00110 bool useEGammaSupercluster) {
00111
00112 mvaEleCut_ = mvaEleCut;
00113 usePFElectrons_ = usePFElectrons;
00114 applyCrackCorrectionsElectrons_ = applyCrackCorrections;
00115 usePFSCEleCalib_ = usePFSCEleCalib;
00116 thePFSCEnergyCalibration_ = thePFSCEnergyCalibration;
00117 useEGElectrons_ = useEGElectrons;
00118 useEGammaSupercluster_ = useEGammaSupercluster;
00119 sumEtEcalIsoForEgammaSC_barrel_ = sumEtEcalIsoForEgammaSC_barrel;
00120 sumEtEcalIsoForEgammaSC_endcap_ = sumEtEcalIsoForEgammaSC_endcap;
00121 coneEcalIsoForEgammaSC_ = coneEcalIsoForEgammaSC;
00122 sumPtTrackIsoForEgammaSC_barrel_ = sumPtTrackIsoForEgammaSC_barrel;
00123 sumPtTrackIsoForEgammaSC_endcap_ = sumPtTrackIsoForEgammaSC_endcap;
00124 coneTrackIsoForEgammaSC_ = coneTrackIsoForEgammaSC;
00125 nTrackIsoForEgammaSC_ = nTrackIsoForEgammaSC;
00126
00127
00128 if(!usePFElectrons_) return;
00129 mvaWeightFileEleID_ = mvaWeightFileEleID;
00130 FILE * fileEleID = fopen(mvaWeightFileEleID_.c_str(), "r");
00131 if (fileEleID) {
00132 fclose(fileEleID);
00133 }
00134 else {
00135 string err = "PFAlgo: cannot open weight file '";
00136 err += mvaWeightFileEleID;
00137 err += "'";
00138 throw invalid_argument( err );
00139 }
00140 pfele_= new PFElectronAlgo(mvaEleCut_,mvaWeightFileEleID_,
00141 thePFSCEnergyCalibration_,
00142 thePFEnergyCalibration,
00143 applyCrackCorrectionsElectrons_,
00144 usePFSCEleCalib_,
00145 useEGElectrons_,
00146 useEGammaSupercluster_,
00147 sumEtEcalIsoForEgammaSC_barrel_,
00148 sumEtEcalIsoForEgammaSC_endcap_,
00149 coneEcalIsoForEgammaSC_,
00150 sumPtTrackIsoForEgammaSC_barrel_,
00151 sumPtTrackIsoForEgammaSC_endcap_,
00152 nTrackIsoForEgammaSC_,
00153 coneTrackIsoForEgammaSC_);
00154 }
00155
00156 void
00157 PFAlgo::setPFPhotonParameters(bool usePFPhotons,
00158 std::string mvaWeightFileConvID,
00159 double mvaConvCut,
00160 bool useReg,
00161 std::string X0_Map,
00162 const boost::shared_ptr<PFEnergyCalibration>& thePFEnergyCalibration,
00163 double sumPtTrackIsoForPhoton,
00164 double sumPtTrackIsoSlopeForPhoton)
00165 {
00166
00167 usePFPhotons_ = usePFPhotons;
00168
00169
00170 reco::Vertex dummy;
00171 if(useVertices_)
00172 {
00173 dummy = primaryVertex_;
00174 }
00175 else {
00176 reco::Vertex::Error e;
00177 e(0, 0) = 0.0015 * 0.0015;
00178 e(1, 1) = 0.0015 * 0.0015;
00179 e(2, 2) = 15. * 15.;
00180 reco::Vertex::Point p(0, 0, 0);
00181 dummy = reco::Vertex(p, e, 0, 0, 0);
00182 }
00183
00184 if(! usePFPhotons_) return;
00185 FILE * filePhotonConvID = fopen(mvaWeightFileConvID.c_str(), "r");
00186 if (filePhotonConvID) {
00187 fclose(filePhotonConvID);
00188 }
00189 else {
00190 string err = "PFAlgo: cannot open weight file '";
00191 err += mvaWeightFileConvID;
00192 err += "'";
00193 throw invalid_argument( err );
00194 }
00195 const reco::Vertex* pv=&dummy;
00196 pfpho_ = new PFPhotonAlgo(mvaWeightFileConvID,
00197 mvaConvCut,
00198 useReg,
00199 X0_Map,
00200 *pv,
00201 thePFEnergyCalibration,
00202 sumPtTrackIsoForPhoton,
00203 sumPtTrackIsoSlopeForPhoton
00204 );
00205 return;
00206 }
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217 void PFAlgo::setPFPhotonRegWeights(
00218 const GBRForest *LCorrForestEB,
00219 const GBRForest *LCorrForestEE,
00220 const GBRForest *GCorrForestBarrel,
00221 const GBRForest *GCorrForestEndcapHr9,
00222 const GBRForest *GCorrForestEndcapLr9, const GBRForest *PFEcalResolution
00223 ){
00224
00225 pfpho_->setGBRForest(LCorrForestEB,LCorrForestEE,
00226 GCorrForestBarrel, GCorrForestEndcapHr9,
00227 GCorrForestEndcapLr9, PFEcalResolution);
00228 }
00229 void
00230 PFAlgo::setPFMuonAndFakeParameters(const edm::ParameterSet& pset)
00231
00232
00233
00234
00235
00236
00237 {
00238
00239
00240
00241
00242 pfmu_ = new PFMuonAlgo();
00243 pfmu_->setParameters(pset);
00244
00245
00246 muonHCAL_= pset.getParameter<std::vector<double> >("muon_HCAL");
00247 muonECAL_= pset.getParameter<std::vector<double> >("muon_ECAL");
00248 muonHO_= pset.getParameter<std::vector<double> >("muon_HO");
00249 assert ( muonHCAL_.size() == 2 && muonECAL_.size() == 2 && muonHO_.size() == 2);
00250 nSigmaTRACK_= pset.getParameter<double>("nsigma_TRACK");
00251 ptError_= pset.getParameter<double>("pt_Error");
00252 factors45_ = pset.getParameter<std::vector<double> >("factors_45");
00253 assert ( factors45_.size() == 2 );
00254
00255 }
00256
00257
00258 void
00259 PFAlgo::setMuonHandle(const edm::Handle<reco::MuonCollection>& muons) {
00260 muonHandle_ = muons;
00261 }
00262
00263
00264 void
00265 PFAlgo::setPostHFCleaningParameters(bool postHFCleaning,
00266 double minHFCleaningPt,
00267 double minSignificance,
00268 double maxSignificance,
00269 double minSignificanceReduction,
00270 double maxDeltaPhiPt,
00271 double minDeltaMet) {
00272 postHFCleaning_ = postHFCleaning;
00273 minHFCleaningPt_ = minHFCleaningPt;
00274 minSignificance_ = minSignificance;
00275 maxSignificance_ = maxSignificance;
00276 minSignificanceReduction_= minSignificanceReduction;
00277 maxDeltaPhiPt_ = maxDeltaPhiPt;
00278 minDeltaMet_ = minDeltaMet;
00279 }
00280
00281 void
00282 PFAlgo::setDisplacedVerticesParameters(bool rejectTracks_Bad,
00283 bool rejectTracks_Step45,
00284 bool usePFNuclearInteractions,
00285 bool usePFConversions,
00286 bool usePFDecays,
00287 double dptRel_DispVtx){
00288
00289 rejectTracks_Bad_ = rejectTracks_Bad;
00290 rejectTracks_Step45_ = rejectTracks_Step45;
00291 usePFNuclearInteractions_ = usePFNuclearInteractions;
00292 usePFConversions_ = usePFConversions;
00293 usePFDecays_ = usePFDecays;
00294 dptRel_DispVtx_ = dptRel_DispVtx;
00295
00296 }
00297
00298
00299 void
00300 PFAlgo::setPFVertexParameters(bool useVertex,
00301 const reco::VertexCollection* primaryVertices) {
00302 useVertices_ = useVertex;
00303
00304
00305 pfmu_->setInputsForCleaning(primaryVertices);
00306
00307
00308
00309 bool primaryVertexFound = false;
00310 int nVtx=primaryVertices->size();
00311 if(usePFPhotons_){
00312 pfpho_->setnPU(nVtx);
00313 }
00314 for (unsigned short i=0 ;i<primaryVertices->size();++i)
00315 {
00316 if(primaryVertices->at(i).isValid()&&(!primaryVertices->at(i).isFake()))
00317 {
00318 primaryVertex_ = primaryVertices->at(i);
00319 primaryVertexFound = true;
00320 break;
00321 }
00322 }
00323
00324 useVertices_ = useVertex && primaryVertexFound;
00325 if(usePFPhotons_) {
00326 if (useVertices_ ){
00327 pfpho_->setPhotonPrimaryVtx(primaryVertex_ );
00328 }
00329 else{
00330 reco::Vertex::Error e;
00331 e(0, 0) = 0.0015 * 0.0015;
00332 e(1, 1) = 0.0015 * 0.0015;
00333 e(2, 2) = 15. * 15.;
00334 reco::Vertex::Point p(0, 0, 0);
00335 reco::Vertex dummy = reco::Vertex(p, e, 0, 0, 0);
00336
00337 pfpho_->setPhotonPrimaryVtx(dummy);
00338 }
00339 }
00340 }
00341
00342
00343 void PFAlgo::reconstructParticles( const reco::PFBlockHandle& blockHandle ) {
00344
00345 blockHandle_ = blockHandle;
00346 reconstructParticles( *blockHandle_ );
00347 }
00348
00349
00350
00351 void PFAlgo::reconstructParticles( const reco::PFBlockCollection& blocks ) {
00352
00353
00354 if(pfCandidates_.get() )
00355 pfCandidates_->clear();
00356 else
00357 pfCandidates_.reset( new reco::PFCandidateCollection );
00358
00359 if(pfElectronCandidates_.get() )
00360 pfElectronCandidates_->clear();
00361 else
00362 pfElectronCandidates_.reset( new reco::PFCandidateCollection);
00363
00364
00365 if( pfPhotonCandidates_.get() )
00366 pfPhotonCandidates_->clear();
00367 else
00368 pfPhotonCandidates_.reset( new reco::PFCandidateCollection);
00369
00370 if(pfCleanedCandidates_.get() )
00371 pfCleanedCandidates_->clear();
00372 else
00373 pfCleanedCandidates_.reset( new reco::PFCandidateCollection );
00374
00375
00376 pfElectronExtra_.clear();
00377 pfPhotonExtra_.clear();
00378
00379 if( debug_ ) {
00380 cout<<"*********************************************************"<<endl;
00381 cout<<"***** Particle flow algorithm *****"<<endl;
00382 cout<<"*********************************************************"<<endl;
00383 }
00384
00385
00386 std::list< reco::PFBlockRef > hcalBlockRefs;
00387 std::list< reco::PFBlockRef > ecalBlockRefs;
00388 std::list< reco::PFBlockRef > hoBlockRefs;
00389 std::list< reco::PFBlockRef > otherBlockRefs;
00390
00391 for( unsigned i=0; i<blocks.size(); ++i ) {
00392
00393 reco::PFBlockRef blockref = createBlockRef( blocks, i);
00394
00395 const reco::PFBlock& block = *blockref;
00396 const edm::OwnVector< reco::PFBlockElement >&
00397 elements = block.elements();
00398
00399 bool singleEcalOrHcal = false;
00400 if( elements.size() == 1 ){
00401 if( elements[0].type() == reco::PFBlockElement::ECAL ){
00402 ecalBlockRefs.push_back( blockref );
00403 singleEcalOrHcal = true;
00404 }
00405 if( elements[0].type() == reco::PFBlockElement::HCAL ){
00406 hcalBlockRefs.push_back( blockref );
00407 singleEcalOrHcal = true;
00408 }
00409 if( elements[0].type() == reco::PFBlockElement::HO ){
00410
00411 hoBlockRefs.push_back( blockref );
00412 singleEcalOrHcal = true;
00413 }
00414 }
00415
00416 if(!singleEcalOrHcal) {
00417 otherBlockRefs.push_back( blockref );
00418 }
00419 }
00420
00421 if( debug_ ){
00422 cout<<"# Ecal blocks: "<<ecalBlockRefs.size()
00423 <<", # Hcal blocks: "<<hcalBlockRefs.size()
00424 <<", # HO blocks: "<<hoBlockRefs.size()
00425 <<", # Other blocks: "<<otherBlockRefs.size()<<endl;
00426 }
00427
00428
00429
00430
00431
00432 unsigned nblcks = 0;
00433 for( IBR io = otherBlockRefs.begin(); io!=otherBlockRefs.end(); ++io) {
00434 if ( debug_ ) std::cout << "Block number " << nblcks++ << std::endl;
00435 processBlock( *io, hcalBlockRefs, ecalBlockRefs );
00436 }
00437
00438 std::list< reco::PFBlockRef > empty;
00439
00440 unsigned hblcks = 0;
00441
00442 for( IBR ih = hcalBlockRefs.begin(); ih!=hcalBlockRefs.end(); ++ih) {
00443 if ( debug_ ) std::cout << "HCAL block number " << hblcks++ << std::endl;
00444 processBlock( *ih, empty, empty );
00445 }
00446
00447 unsigned eblcks = 0;
00448
00449 for( IBR ie = ecalBlockRefs.begin(); ie!=ecalBlockRefs.end(); ++ie) {
00450 if ( debug_ ) std::cout << "ECAL block number " << eblcks++ << std::endl;
00451 processBlock( *ie, empty, empty );
00452 }
00453
00454
00455 postCleaning();
00456
00457
00458 pfmu_->postClean(pfCandidates_.get());
00459
00460
00461 if( muonHandle_.isValid())
00462 pfmu_->addMissingMuons(muonHandle_,pfCandidates_.get());
00463 }
00464
00465
00466 void PFAlgo::processBlock( const reco::PFBlockRef& blockref,
00467 std::list<reco::PFBlockRef>& hcalBlockRefs,
00468 std::list<reco::PFBlockRef>& ecalBlockRefs ) {
00469
00470
00471 assert(!blockref.isNull() );
00472 const reco::PFBlock& block = *blockref;
00473
00474 typedef std::multimap<double, unsigned>::iterator IE;
00475 typedef std::multimap<double, std::pair<unsigned,math::XYZVector> >::iterator IS;
00476 typedef std::multimap<double, std::pair<unsigned,bool> >::iterator IT;
00477 typedef std::multimap< unsigned, std::pair<double, unsigned> >::iterator II;
00478
00479 if(debug_) {
00480 cout<<"#########################################################"<<endl;
00481 cout<<"##### Process Block: #####"<<endl;
00482 cout<<"#########################################################"<<endl;
00483 cout<<block<<endl;
00484 }
00485
00486
00487 const edm::OwnVector< reco::PFBlockElement >& elements = block.elements();
00488
00489 PFBlock::LinkData linkData = block.linkData();
00490
00491
00492 vector<bool> active( elements.size(), true );
00493
00494
00495
00496
00497 std::vector<reco::PFCandidate> tempElectronCandidates;
00498 tempElectronCandidates.clear();
00499 if (usePFElectrons_) {
00500 if (pfele_->isElectronValidCandidate(blockref,active, primaryVertex_ )){
00501
00502 const std::vector<reco::PFCandidate> PFElectCandidates_(pfele_->getElectronCandidates());
00503 for ( std::vector<reco::PFCandidate>::const_iterator ec=PFElectCandidates_.begin(); ec != PFElectCandidates_.end(); ++ec )tempElectronCandidates.push_back(*ec);
00504
00505
00506
00507
00508 }
00509
00510
00511 pfElectronCandidates_->insert(pfElectronCandidates_->end(),
00512 pfele_->getAllElectronCandidates().begin(),
00513 pfele_->getAllElectronCandidates().end());
00514
00515 pfElectronExtra_.insert(pfElectronExtra_.end(),
00516 pfele_->getElectronExtra().begin(),
00517 pfele_->getElectronExtra().end());
00518
00519 }
00520 if( usePFPhotons_ ) {
00521
00522 if(debug_)
00523 cout<<endl<<"--------------- entering PFPhotonAlgo ----------------"<<endl;
00524 vector<PFCandidatePhotonExtra> pfPhotonExtraCand;
00525 if ( pfpho_->isPhotonValidCandidate(blockref,
00526 active,
00527 pfPhotonCandidates_,
00528 pfPhotonExtraCand,
00529 tempElectronCandidates
00530
00531 ) ) {
00532 if(debug_)
00533 std::cout<< " In this PFBlock we found "<<pfPhotonCandidates_->size()<<" Photon Candidates."<<std::endl;
00534
00535
00536
00537
00538
00539 unsigned int extracand =0;
00540 PFCandidateCollection::const_iterator cand = pfPhotonCandidates_->begin();
00541 for( ; cand != pfPhotonCandidates_->end(); ++cand, ++extracand) {
00542 pfCandidates_->push_back(*cand);
00543 pfPhotonExtra_.push_back(pfPhotonExtraCand[extracand]);
00544 }
00545
00546 }
00547 pfPhotonExtraCand.clear();
00548 pfPhotonCandidates_->clear();
00549 }
00550
00551
00552 if (usePFConversions_ )
00553 {
00554 for(unsigned iEle=0; iEle<elements.size(); iEle++) {
00555 PFBlockElement::Type type = elements[iEle].type();
00556 if(type==PFBlockElement::TRACK)
00557 {
00558 if(elements[iEle].trackRef()->algo() == 12)
00559 active[iEle]=false;
00560 if(elements[iEle].trackRef()->quality(reco::TrackBase::highPurity))continue;
00561 const reco::PFBlockElementTrack * trackRef = dynamic_cast<const reco::PFBlockElementTrack*>((&elements[iEle]));
00562 if(!(trackRef->trackType(reco::PFBlockElement::T_FROM_GAMMACONV)))continue;
00563 if(elements[iEle].convRef().isNonnull())active[iEle]=false;
00564 }
00565 }
00566 }
00567 for ( std::vector<reco::PFCandidate>::const_iterator ec=tempElectronCandidates.begin(); ec != tempElectronCandidates.end(); ++ec ){
00568 pfCandidates_->push_back(*ec);
00569 }
00570 tempElectronCandidates.clear();
00571
00572
00573 if(debug_)
00574 cout<<endl<<"--------------- loop 1 ------------------"<<endl;
00575
00576
00577
00578
00579
00580
00581
00582
00583
00584
00585
00586
00587
00588
00589
00590
00591
00592
00593
00594
00595
00596
00597
00598
00599
00600
00601 vector<unsigned> hcalIs;
00602 vector<unsigned> hoIs;
00603 vector<unsigned> ecalIs;
00604 vector<unsigned> trackIs;
00605 vector<unsigned> ps1Is;
00606 vector<unsigned> ps2Is;
00607
00608 vector<unsigned> hfEmIs;
00609 vector<unsigned> hfHadIs;
00610
00611
00612 for(unsigned iEle=0; iEle<elements.size(); iEle++) {
00613 PFBlockElement::Type type = elements[iEle].type();
00614
00615 if(debug_ && type != PFBlockElement::BREM ) cout<<endl<<elements[iEle];
00616
00617 switch( type ) {
00618 case PFBlockElement::TRACK:
00619 if ( active[iEle] ) {
00620 trackIs.push_back( iEle );
00621 if(debug_) cout<<"TRACK, stored index, continue"<<endl;
00622 }
00623 break;
00624 case PFBlockElement::ECAL:
00625 if ( active[iEle] ) {
00626 ecalIs.push_back( iEle );
00627 if(debug_) cout<<"ECAL, stored index, continue"<<endl;
00628 }
00629 continue;
00630 case PFBlockElement::HCAL:
00631 if ( active[iEle] ) {
00632 hcalIs.push_back( iEle );
00633 if(debug_) cout<<"HCAL, stored index, continue"<<endl;
00634 }
00635 continue;
00636 case PFBlockElement::HO:
00637 if (useHO_) {
00638 if ( active[iEle] ) {
00639 hoIs.push_back( iEle );
00640 if(debug_) cout<<"HO, stored index, continue"<<endl;
00641 }
00642 }
00643 continue;
00644 case PFBlockElement::HFEM:
00645 if ( active[iEle] ) {
00646 hfEmIs.push_back( iEle );
00647 if(debug_) cout<<"HFEM, stored index, continue"<<endl;
00648 }
00649 continue;
00650 case PFBlockElement::HFHAD:
00651 if ( active[iEle] ) {
00652 hfHadIs.push_back( iEle );
00653 if(debug_) cout<<"HFHAD, stored index, continue"<<endl;
00654 }
00655 continue;
00656 default:
00657 continue;
00658 }
00659
00660
00661 unsigned iTrack = iEle;
00662 reco::MuonRef muonRef = elements[iTrack].muonRef();
00663
00664
00665
00666
00667 if (active[iTrack] && isFromSecInt(elements[iEle], "primary")){
00668 bool isPrimaryTrack = elements[iEle].displacedVertexRef(PFBlockElement::T_TO_DISP)->displacedVertexRef()->isTherePrimaryTracks();
00669 if (isPrimaryTrack) {
00670 if (debug_) cout << "Primary Track reconstructed alone" << endl;
00671
00672 unsigned tmpi = reconstructTrack(elements[iEle]);
00673 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iEle );
00674 active[iTrack] = false;
00675 }
00676 }
00677
00678
00679 if(debug_) {
00680 if ( !active[iTrack] )
00681 cout << "Already used by electrons, muons, conversions" << endl;
00682 }
00683
00684
00685
00686 if ( ! active[iTrack] ) continue;
00687
00688 reco::TrackRef trackRef = elements[iTrack].trackRef();
00689 assert( !trackRef.isNull() );
00690
00691 if (debug_ ) {
00692 cout <<"PFAlgo:processBlock "<<" "<<trackIs.size()<<" "<<ecalIs.size()<<" "<<hcalIs.size()<<" "<<hoIs.size()<<endl;
00693 }
00694
00695
00696
00697
00698
00699 std::multimap<double, unsigned> ecalElems;
00700 block.associatedElements( iTrack, linkData,
00701 ecalElems ,
00702 reco::PFBlockElement::ECAL,
00703 reco::PFBlock::LINKTEST_ALL );
00704
00705 std::multimap<double, unsigned> hcalElems;
00706 block.associatedElements( iTrack, linkData,
00707 hcalElems,
00708 reco::PFBlockElement::HCAL,
00709 reco::PFBlock::LINKTEST_ALL );
00710
00711
00712
00713
00714 if ( hcalElems.empty() && !ecalElems.empty() ) {
00715
00716 unsigned ntt = 1;
00717 unsigned index = ecalElems.begin()->second;
00718 std::multimap<double, unsigned> sortedTracks;
00719 block.associatedElements( index, linkData,
00720 sortedTracks,
00721 reco::PFBlockElement::TRACK,
00722 reco::PFBlock::LINKTEST_ALL );
00723
00724
00725
00726 for(IE ie = sortedTracks.begin(); ie != sortedTracks.end(); ++ie ) {
00727 unsigned jTrack = ie->second;
00728
00729
00730 if ( !active[jTrack] ) continue;
00731
00732
00733
00734 if ( jTrack == iTrack ) continue;
00735
00736
00737
00738
00739 std::multimap<double, unsigned> sortedECAL;
00740 block.associatedElements( jTrack, linkData,
00741 sortedECAL,
00742 reco::PFBlockElement::ECAL,
00743 reco::PFBlock::LINKTEST_ALL );
00744 if ( sortedECAL.begin()->second != index ) continue;
00745
00746
00747
00748
00749 std::multimap<double, unsigned> sortedHCAL;
00750 block.associatedElements( jTrack, linkData,
00751 sortedHCAL,
00752 reco::PFBlockElement::HCAL,
00753 reco::PFBlock::LINKTEST_ALL );
00754 if ( !sortedHCAL.size() ) continue;
00755
00756 ntt++;
00757
00758
00759 block.setLink( iTrack,
00760 sortedHCAL.begin()->second,
00761 sortedECAL.begin()->first,
00762 linkData,
00763 PFBlock::LINKTEST_RECHIT );
00764
00765 }
00766
00767
00768 block.associatedElements( iTrack, linkData,
00769 hcalElems,
00770 reco::PFBlockElement::HCAL,
00771 reco::PFBlock::LINKTEST_ALL );
00772
00773 if ( debug_ && hcalElems.size() )
00774 std::cout << "Track linked back to HCAL due to ECAL sharing with other tracks" << std::endl;
00775 }
00776
00777
00778
00779
00780
00781
00782
00783
00784 std::multimap<double,unsigned> gsfElems;
00785 if (usePFElectrons_ == false) {
00786 block.associatedElements( iTrack, linkData,
00787 gsfElems ,
00788 reco::PFBlockElement::GSF );
00789 }
00790
00791
00792 if(hcalElems.empty() && debug_) {
00793 cout<<"no hcal element connected to track "<<iTrack<<endl;
00794 }
00795
00796
00797
00798
00799 bool hcalFound = false;
00800
00801 if(debug_)
00802 cout<<"now looping on elements associated to the track"<<endl;
00803
00804
00805
00806 for(IE ie = ecalElems.begin(); ie != ecalElems.end(); ++ie ) {
00807
00808 unsigned index = ie->second;
00809
00810 PFBlockElement::Type type = elements[index].type();
00811 if(debug_) {
00812 double dist = ie->first;
00813 cout<<"\telement "<<elements[index]<<" linked with distance = "<< dist <<endl;
00814 if ( ! active[index] ) cout << "This ECAL is already used - skip it" << endl;
00815 }
00816 assert( type == PFBlockElement::ECAL );
00817
00818
00819 if ( ! active[index] ) continue;
00820
00821
00822
00823
00824
00825 if( !hcalElems.empty() && debug_)
00826 cout<<"\t\tat least one hcal element connected to the track."
00827 <<" Sparing Ecal cluster for the hcal loop"<<endl;
00828
00829 }
00830
00831
00832
00833
00834
00835 if( hcalElems.empty() ) {
00836
00837 if ( debug_ )
00838 std::cout << "Now deals with tracks linked to no HCAL clusters" << std::endl;
00839
00840
00841
00842
00843 double trackMomentum = elements[iTrack].trackRef()->p();
00844 if ( debug_ )
00845 std::cout << elements[iTrack] << std::endl;
00846
00847
00848
00849
00850 bool thisIsAMuon = PFMuonAlgo::isMuon(elements[iTrack]);
00851 if ( thisIsAMuon ) trackMomentum = 0.;
00852
00853
00854
00855 bool rejectFake = false;
00856 if ( !thisIsAMuon && elements[iTrack].trackRef()->ptError() > ptError_ ) {
00857
00858 double deficit = trackMomentum;
00859 double resol = neutralHadronEnergyResolution(trackMomentum,
00860 elements[iTrack].trackRef()->eta());
00861 resol *= trackMomentum;
00862
00863 if ( !ecalElems.empty() ) {
00864 unsigned thisEcal = ecalElems.begin()->second;
00865 reco::PFClusterRef clusterRef = elements[thisEcal].clusterRef();
00866 deficit -= clusterRef->energy();
00867 resol = neutralHadronEnergyResolution(trackMomentum,
00868 clusterRef->positionREP().Eta());
00869 resol *= trackMomentum;
00870 }
00871
00872 bool isPrimary = isFromSecInt(elements[iTrack], "primary");
00873
00874 if ( deficit > nSigmaTRACK_*resol && !isPrimary) {
00875 rejectFake = true;
00876 active[iTrack] = false;
00877 if ( debug_ )
00878 std::cout << elements[iTrack] << std::endl
00879 << "is probably a fake (1) --> lock the track"
00880 << std::endl;
00881 }
00882 }
00883
00884 if ( rejectFake ) continue;
00885
00886
00887
00888
00889 std::vector<unsigned> tmpi;
00890 std::vector<unsigned> kTrack;
00891
00892
00893 double Dpt = trackRef->ptError();
00894
00895 if ( rejectTracks_Step45_ && ecalElems.empty() &&
00896 trackMomentum > 30. && Dpt > 0.5 &&
00897 ( trackRef->algo() == TrackBase::iter4 ||
00898 trackRef->algo() == TrackBase::iter5 ||
00899 trackRef->algo() == TrackBase::iter6 ) ) {
00900
00901
00902 double dptRel = Dpt/trackRef->pt()*100;
00903 bool isPrimaryOrSecondary = isFromSecInt(elements[iTrack], "all");
00904
00905 if ( isPrimaryOrSecondary && dptRel < dptRel_DispVtx_) continue;
00906 unsigned nHits = elements[iTrack].trackRef()->hitPattern().trackerLayersWithMeasurement();
00907 unsigned int NLostHit = trackRef->hitPattern().trackerLayersWithoutMeasurement();
00908
00909 if ( debug_ )
00910 std::cout << "A track (algo = " << trackRef->algo() << ") with momentum " << trackMomentum
00911 << " / " << elements[iTrack].trackRef()->pt() << " +/- " << Dpt
00912 << " / " << elements[iTrack].trackRef()->eta()
00913 << " without any link to ECAL/HCAL and with " << nHits << " (" << NLostHit
00914 << ") hits (lost hits) has been cleaned" << std::endl;
00915 active[iTrack] = false;
00916 continue;
00917 }
00918
00919
00920 tmpi.push_back(reconstructTrack( elements[iTrack]));
00921
00922 kTrack.push_back(iTrack);
00923 active[iTrack] = false;
00924
00925
00926 if ( ecalElems.empty() ) {
00927 (*pfCandidates_)[tmpi[0]].setEcalEnergy( 0., 0. );
00928 (*pfCandidates_)[tmpi[0]].setHcalEnergy( 0., 0. );
00929 (*pfCandidates_)[tmpi[0]].setHoEnergy( 0., 0. );
00930 (*pfCandidates_)[tmpi[0]].setPs1Energy( 0 );
00931 (*pfCandidates_)[tmpi[0]].setPs2Energy( 0 );
00932 (*pfCandidates_)[tmpi[0]].addElementInBlock( blockref, kTrack[0] );
00933 continue;
00934 }
00935
00936
00937 unsigned thisEcal = ecalElems.begin()->second;
00938 reco::PFClusterRef clusterRef = elements[thisEcal].clusterRef();
00939 if ( debug_ ) std::cout << " is associated to " << elements[thisEcal] << std::endl;
00940
00941
00942
00943 if ( thisIsAMuon ) {
00944 (*pfCandidates_)[tmpi[0]].setEcalEnergy( clusterRef->energy(),
00945 std::min(clusterRef->energy(), muonECAL_[0]) );
00946 (*pfCandidates_)[tmpi[0]].setHcalEnergy( 0., 0. );
00947 (*pfCandidates_)[tmpi[0]].setHoEnergy( 0., 0. );
00948 (*pfCandidates_)[tmpi[0]].setPs1Energy( 0 );
00949 (*pfCandidates_)[tmpi[0]].setPs2Energy( 0 );
00950 (*pfCandidates_)[tmpi[0]].addElementInBlock( blockref, kTrack[0] );
00951 }
00952
00953 double slopeEcal = 1.;
00954 bool connectedToEcal = false;
00955 unsigned iEcal = 99999;
00956 double calibEcal = 0.;
00957 double calibHcal = 0.;
00958 double totalEcal = thisIsAMuon ? -muonECAL_[0] : 0.;
00959
00960
00961 std::multimap<double, unsigned> sortedTracks;
00962 block.associatedElements( thisEcal, linkData,
00963 sortedTracks,
00964 reco::PFBlockElement::TRACK,
00965 reco::PFBlock::LINKTEST_ALL );
00966
00967 for(IE ie = sortedTracks.begin(); ie != sortedTracks.end(); ++ie ) {
00968 unsigned jTrack = ie->second;
00969
00970
00971 if ( !active[jTrack] ) continue;
00972
00973
00974 if ( jTrack == iTrack ) continue;
00975
00976
00977
00978
00979 std::multimap<double, unsigned> sortedECAL;
00980 block.associatedElements( jTrack, linkData,
00981 sortedECAL,
00982 reco::PFBlockElement::ECAL,
00983 reco::PFBlock::LINKTEST_ALL );
00984 if ( sortedECAL.begin()->second != thisEcal ) continue;
00985
00986
00987 bool thatIsAMuon = PFMuonAlgo::isMuon(elements[jTrack]);
00988
00989
00990
00991 bool rejectFake = false;
00992 reco::TrackRef trackRef = elements[jTrack].trackRef();
00993 if ( !thatIsAMuon && trackRef->ptError() > ptError_) {
00994 double deficit = trackMomentum + trackRef->p() - clusterRef->energy();
00995 double resol = nSigmaTRACK_*neutralHadronEnergyResolution(trackMomentum+trackRef->p(),
00996 clusterRef->positionREP().Eta());
00997 resol *= (trackMomentum+trackRef->p());
00998 if ( deficit > nSigmaTRACK_*resol ) {
00999 rejectFake = true;
01000 kTrack.push_back(jTrack);
01001 active[jTrack] = false;
01002 if ( debug_ )
01003 std::cout << elements[jTrack] << std::endl
01004 << "is probably a fake (2) --> lock the track"
01005 << std::endl;
01006 }
01007 }
01008 if ( rejectFake ) continue;
01009
01010
01011
01012
01013
01014 if ( !thatIsAMuon ) {
01015 if ( debug_ )
01016 std::cout << "Track momentum increased from " << trackMomentum << " GeV ";
01017 trackMomentum += trackRef->p();
01018 if ( debug_ ) {
01019 std::cout << "to " << trackMomentum << " GeV." << std::endl;
01020 std::cout << "with " << elements[jTrack] << std::endl;
01021 }
01022 } else {
01023 totalEcal -= muonECAL_[0];
01024 totalEcal = std::max(totalEcal, 0.);
01025 }
01026
01027
01028
01029 tmpi.push_back(reconstructTrack( elements[jTrack] ));
01030
01031
01032 kTrack.push_back(jTrack);
01033 active[jTrack] = false;
01034
01035 if ( thatIsAMuon ) {
01036 (*pfCandidates_)[tmpi.back()].setEcalEnergy(clusterRef->energy(),
01037 std::min(clusterRef->energy(),muonECAL_[0]));
01038 (*pfCandidates_)[tmpi.back()].setHcalEnergy( 0., 0. );
01039 (*pfCandidates_)[tmpi.back()].setHoEnergy( 0., 0. );
01040 (*pfCandidates_)[tmpi.back()].setPs1Energy( 0 );
01041 (*pfCandidates_)[tmpi.back()].setPs2Energy( 0 );
01042 (*pfCandidates_)[tmpi.back()].addElementInBlock( blockref, kTrack.back() );
01043 }
01044 }
01045
01046
01047 if ( debug_ ) std::cout << "Loop over all associated ECAL clusters" << std::endl;
01048
01049 for(IE ie = ecalElems.begin(); ie != ecalElems.end(); ++ie ) {
01050
01051 unsigned index = ie->second;
01052 PFBlockElement::Type type = elements[index].type();
01053 assert( type == PFBlockElement::ECAL );
01054 if ( debug_ ) std::cout << elements[index] << std::endl;
01055
01056
01057 if ( debug_ && ! active[index] ) std::cout << "is not active - ignore " << std::endl;
01058 if ( ! active[index] ) continue;
01059
01060
01061
01062 block.associatedElements( index, linkData,
01063 sortedTracks,
01064 reco::PFBlockElement::TRACK,
01065 reco::PFBlock::LINKTEST_ALL );
01066 bool skip = true;
01067 for (unsigned ic=0; ic<kTrack.size();++ic) {
01068 if ( sortedTracks.begin()->second == kTrack[ic] ) {
01069 skip = false;
01070 break;
01071 }
01072 }
01073 if ( debug_ && skip ) std::cout << "is closer to another track - ignore " << std::endl;
01074 if ( skip ) continue;
01075
01076
01077
01078 reco::PFClusterRef clusterRef = elements[index].clusterRef();
01079 assert( !clusterRef.isNull() );
01080
01081 if ( debug_ ) {
01082 double dist = ie->first;
01083 std::cout <<"Ecal cluster with raw energy = " << clusterRef->energy()
01084 <<" linked with distance = " << dist << std::endl;
01085 }
01086
01087
01088
01089
01090
01091
01092
01093
01094
01095
01096
01097 vector<double> ps1Ene(1,static_cast<double>(0.));
01098 associatePSClusters(index, reco::PFBlockElement::PS1, block, elements, linkData, active, ps1Ene);
01099 vector<double> ps2Ene(1,static_cast<double>(0.));
01100 associatePSClusters(index, reco::PFBlockElement::PS2, block, elements, linkData, active, ps2Ene);
01101
01102
01103 bool crackCorrection = false;
01104 double ecalEnergy = calibration_->energyEm(*clusterRef,ps1Ene,ps2Ene,crackCorrection);
01105 if ( debug_ )
01106 std::cout << "Corrected ECAL(+PS) energy = " << ecalEnergy << std::endl;
01107
01108
01109
01110 totalEcal += ecalEnergy;
01111 double previousCalibEcal = calibEcal;
01112 double previousSlopeEcal = slopeEcal;
01113 calibEcal = std::max(totalEcal,0.);
01114 calibHcal = 0.;
01115 calibration_->energyEmHad(trackMomentum,calibEcal,calibHcal,
01116 clusterRef->positionREP().Eta(),
01117 clusterRef->positionREP().Phi());
01118 if ( totalEcal > 0.) slopeEcal = calibEcal/totalEcal;
01119
01120 if ( debug_ )
01121 std::cout << "The total calibrated energy so far amounts to = " << calibEcal << std::endl;
01122
01123
01124
01125 if ( connectedToEcal && calibEcal - trackMomentum >= 0. ) {
01126
01127
01128 calibEcal = previousCalibEcal;
01129 slopeEcal = previousSlopeEcal;
01130 totalEcal = calibEcal/slopeEcal;
01131
01132
01133
01134 active[index] = false;
01135
01136
01137 std::multimap<double, unsigned> assTracks;
01138 block.associatedElements( index, linkData,
01139 assTracks,
01140 reco::PFBlockElement::TRACK,
01141 reco::PFBlock::LINKTEST_ALL );
01142
01143
01144 unsigned tmpe = reconstructCluster( *clusterRef, ecalEnergy );
01145 (*pfCandidates_)[tmpe].setEcalEnergy( clusterRef->energy(), ecalEnergy );
01146 (*pfCandidates_)[tmpe].setHcalEnergy( 0., 0. );
01147 (*pfCandidates_)[tmpe].setHoEnergy( 0., 0. );
01148 (*pfCandidates_)[tmpe].setPs1Energy( ps1Ene[0] );
01149 (*pfCandidates_)[tmpe].setPs2Energy( ps2Ene[0] );
01150 (*pfCandidates_)[tmpe].addElementInBlock( blockref, index );
01151
01152 if(assTracks.size()) {
01153 (*pfCandidates_)[tmpe].addElementInBlock( blockref, assTracks.begin()->second );
01154
01155
01156 const math::XYZPointF& chargedPosition =
01157 dynamic_cast<const reco::PFBlockElementTrack*>(&elements[assTracks.begin()->second])->positionAtECALEntrance();
01158 (*pfCandidates_)[tmpe].setPositionAtECALEntrance(chargedPosition);
01159 }
01160 break;
01161 }
01162
01163
01164 connectedToEcal = true;
01165 iEcal = index;
01166 active[index] = false;
01167 for (unsigned ic=0; ic<tmpi.size();++ic)
01168 (*pfCandidates_)[tmpi[ic]].addElementInBlock( blockref, iEcal );
01169
01170
01171 }
01172
01173 bool bNeutralProduced = false;
01174
01175
01176 if( connectedToEcal ) {
01177 reco::PFClusterRef pivotalRef = elements[iEcal].clusterRef();
01178
01179 double neutralEnergy = calibEcal - trackMomentum;
01180
01181
01182
01183
01184
01185
01186
01187
01188
01189
01190
01191
01192
01193
01194
01195
01196
01197
01198
01199
01200
01201
01202
01203
01204
01205
01206
01207
01208
01209
01210
01211
01212
01213
01214
01215 double resol = neutralHadronEnergyResolution(trackMomentum,pivotalRef->positionREP().Eta());
01216 resol *= trackMomentum;
01217 if ( neutralEnergy > std::max(0.5,nSigmaECAL_*resol) ) {
01218 neutralEnergy /= slopeEcal;
01219 unsigned tmpj = reconstructCluster( *pivotalRef, neutralEnergy );
01220 (*pfCandidates_)[tmpj].setEcalEnergy( pivotalRef->energy(), neutralEnergy );
01221 (*pfCandidates_)[tmpj].setHcalEnergy( 0., 0. );
01222 (*pfCandidates_)[tmpj].setHoEnergy( 0., 0. );
01223 (*pfCandidates_)[tmpj].setPs1Energy( 0. );
01224 (*pfCandidates_)[tmpj].setPs2Energy( 0. );
01225 (*pfCandidates_)[tmpj].addElementInBlock(blockref, iEcal);
01226 bNeutralProduced = true;
01227 for (unsigned ic=0; ic<kTrack.size();++ic)
01228 (*pfCandidates_)[tmpj].addElementInBlock( blockref, kTrack[ic] );
01229 }
01230
01231
01232 for (unsigned ic=0; ic<tmpi.size();++ic) {
01233
01234
01235 if ( (*pfCandidates_)[tmpi[ic]].particleId() == reco::PFCandidate::mu ) continue;
01236
01237 double fraction = (*pfCandidates_)[tmpi[ic]].trackRef()->p()/trackMomentum;
01238 double ecalCal = bNeutralProduced ?
01239 (calibEcal-neutralEnergy*slopeEcal)*fraction : calibEcal*fraction;
01240 double ecalRaw = totalEcal*fraction;
01241
01242 if (debug_) cout << "The fraction after photon supression is " << fraction << " calibrated ecal = " << ecalCal << endl;
01243
01244 (*pfCandidates_)[tmpi[ic]].setEcalEnergy( ecalRaw, ecalCal );
01245 (*pfCandidates_)[tmpi[ic]].setHcalEnergy( 0., 0. );
01246 (*pfCandidates_)[tmpi[ic]].setHoEnergy( 0., 0. );
01247 (*pfCandidates_)[tmpi[ic]].setPs1Energy( 0 );
01248 (*pfCandidates_)[tmpi[ic]].setPs2Energy( 0 );
01249 (*pfCandidates_)[tmpi[ic]].addElementInBlock( blockref, kTrack[ic] );
01250 }
01251
01252 }
01253
01254
01255 for (unsigned ic=0; ic<tmpi.size();++ic) {
01256 const PFCandidate& pfc = (*pfCandidates_)[tmpi[ic]];
01257 const PFCandidate::ElementsInBlocks& eleInBlocks = pfc.elementsInBlocks();
01258 if ( eleInBlocks.size() == 0 ) {
01259 if ( debug_ )std::cout << "Single track / Fill element in block! " << std::endl;
01260 (*pfCandidates_)[tmpi[ic]].addElementInBlock( blockref, kTrack[ic] );
01261 }
01262 }
01263
01264 }
01265
01266
01267
01268
01269 for(IE ie = hcalElems.begin(); ie != hcalElems.end(); ++ie ) {
01270
01271 unsigned index = ie->second;
01272
01273 PFBlockElement::Type type = elements[index].type();
01274
01275 if(debug_) {
01276 double dist = block.dist(iTrack,index,linkData,reco::PFBlock::LINKTEST_ALL);
01277 cout<<"\telement "<<elements[index]<<" linked with distance "<< dist <<endl;
01278 }
01279
01280 assert( type == PFBlockElement::HCAL );
01281
01282
01283
01284 if( !hcalFound ) {
01285 if(debug_)
01286 cout<<"\t\tclosest hcal cluster, doing nothing"<<endl;
01287
01288 hcalFound = true;
01289
01290
01291
01292 }
01293 else {
01294
01295 if(debug_)
01296 cout<<"\t\tsecondary hcal cluster. unlinking"<<endl;
01297 block.setLink( iTrack, index, -1., linkData,
01298 PFBlock::LINKTEST_RECHIT );
01299 }
01300 }
01301 }
01302
01303
01304
01305
01306 if( !(hfEmIs.empty() && hfHadIs.empty() ) ) {
01307
01308
01309 assert( hfEmIs.size() + hfHadIs.size() == elements.size() );
01310
01311 if( elements.size() == 1 ) {
01312
01313 reco::PFClusterRef clusterRef = elements[0].clusterRef();
01314 double energyHF = 0.;
01315 double uncalibratedenergyHF = 0.;
01316 unsigned tmpi = 0;
01317 switch( clusterRef->layer() ) {
01318 case PFLayer::HF_EM:
01319
01320 energyHF = clusterRef->energy();
01321 uncalibratedenergyHF = energyHF;
01322 if(thepfEnergyCalibrationHF_->getcalibHF_use()==true){
01323 energyHF = thepfEnergyCalibrationHF_->energyEm(uncalibratedenergyHF,
01324 clusterRef->positionREP().Eta(),
01325 clusterRef->positionREP().Phi());
01326 }
01327 tmpi = reconstructCluster( *clusterRef, energyHF );
01328 (*pfCandidates_)[tmpi].setEcalEnergy( uncalibratedenergyHF, energyHF );
01329 (*pfCandidates_)[tmpi].setHcalEnergy( 0., 0.);
01330 (*pfCandidates_)[tmpi].setHoEnergy( 0., 0.);
01331 (*pfCandidates_)[tmpi].setPs1Energy( 0. );
01332 (*pfCandidates_)[tmpi].setPs2Energy( 0. );
01333 (*pfCandidates_)[tmpi].addElementInBlock( blockref, hfEmIs[0] );
01334
01335 break;
01336 case PFLayer::HF_HAD:
01337
01338 energyHF = clusterRef->energy();
01339 uncalibratedenergyHF = energyHF;
01340 if(thepfEnergyCalibrationHF_->getcalibHF_use()==true){
01341 energyHF = thepfEnergyCalibrationHF_->energyHad(uncalibratedenergyHF,
01342 clusterRef->positionREP().Eta(),
01343 clusterRef->positionREP().Phi());
01344 }
01345 tmpi = reconstructCluster( *clusterRef, energyHF );
01346 (*pfCandidates_)[tmpi].setHcalEnergy( uncalibratedenergyHF, energyHF );
01347 (*pfCandidates_)[tmpi].setEcalEnergy( 0., 0.);
01348 (*pfCandidates_)[tmpi].setHoEnergy( 0., 0.);
01349 (*pfCandidates_)[tmpi].setPs1Energy( 0. );
01350 (*pfCandidates_)[tmpi].setPs2Energy( 0. );
01351 (*pfCandidates_)[tmpi].addElementInBlock( blockref, hfHadIs[0] );
01352
01353 break;
01354 default:
01355 assert(0);
01356 }
01357 }
01358 else if( elements.size() == 2 ) {
01359
01360 reco::PFClusterRef c0 = elements[0].clusterRef();
01361 reco::PFClusterRef c1 = elements[1].clusterRef();
01362
01363 reco::PFClusterRef cem = c1;
01364 reco::PFClusterRef chad = c0;
01365
01366 if( c0->layer()== PFLayer::HF_EM ) {
01367 if ( c1->layer()== PFLayer::HF_HAD ) {
01368 cem = c0; chad = c1;
01369 } else {
01370 assert(0);
01371 }
01372
01373 double energyHfEm = cem->energy();
01374 double energyHfHad = chad->energy();
01375 double uncalibratedenergyHFEm = energyHfEm;
01376 double uncalibratedenergyHFHad = energyHfHad;
01377 if(thepfEnergyCalibrationHF_->getcalibHF_use()==true){
01378
01379 energyHfEm = thepfEnergyCalibrationHF_->energyEmHad(uncalibratedenergyHFEm,
01380 0.0,
01381 c0->positionREP().Eta(),
01382 c0->positionREP().Phi());
01383 energyHfHad = thepfEnergyCalibrationHF_->energyEmHad(0.0,
01384 uncalibratedenergyHFHad,
01385 c1->positionREP().Eta(),
01386 c1->positionREP().Phi());
01387 }
01388 unsigned tmpi = reconstructCluster( *chad, energyHfEm+energyHfHad );
01389 (*pfCandidates_)[tmpi].setEcalEnergy( uncalibratedenergyHFEm, energyHfEm );
01390 (*pfCandidates_)[tmpi].setHcalEnergy( uncalibratedenergyHFHad, energyHfHad);
01391 (*pfCandidates_)[tmpi].setHoEnergy( 0., 0.);
01392 (*pfCandidates_)[tmpi].setPs1Energy( 0. );
01393 (*pfCandidates_)[tmpi].setPs2Energy( 0. );
01394 (*pfCandidates_)[tmpi].addElementInBlock( blockref, hfEmIs[0] );
01395 (*pfCandidates_)[tmpi].addElementInBlock( blockref, hfHadIs[0] );
01396
01397
01398 }
01399 else {
01400 cerr<<"Warning: 2 elements, but not 1 HFEM and 1 HFHAD"<<endl;
01401 cerr<<block<<endl;
01402
01403
01404 }
01405 }
01406 else {
01407
01408
01409 cerr<<"Warning: HF, but n elem different from 1 or 2"<<endl;
01410 cerr<<block<<endl;
01411
01412
01413 }
01414 }
01415
01416
01417
01418 if(debug_) {
01419 cout<<endl;
01420 cout<<endl<<"--------------- loop hcal ---------------------"<<endl;
01421 }
01422
01423
01424
01425
01426
01427
01428
01429 for(unsigned i=0; i<hcalIs.size(); i++) {
01430
01431 unsigned iHcal= hcalIs[i];
01432 PFBlockElement::Type type = elements[iHcal].type();
01433
01434 assert( type == PFBlockElement::HCAL );
01435
01436 if(debug_) cout<<endl<<elements[iHcal]<<endl;
01437
01438
01439
01440
01441
01442 std::multimap<double, unsigned> sortedTracks;
01443 block.associatedElements( iHcal, linkData,
01444 sortedTracks,
01445 reco::PFBlockElement::TRACK,
01446 reco::PFBlock::LINKTEST_ALL );
01447
01448 std::multimap< unsigned, std::pair<double, unsigned> > associatedEcals;
01449
01450 std::map< unsigned, std::pair<double, double> > associatedPSs;
01451
01452 std::multimap<double, std::pair<unsigned,bool> > associatedTracks;
01453
01454
01455 std::multimap<double,std::pair<unsigned,math::XYZVector> > ecalSatellites;
01456 std::pair<unsigned,math::XYZVector> fakeSatellite = make_pair(iHcal,math::XYZVector(0.,0.,0.));
01457 ecalSatellites.insert( make_pair(-1., fakeSatellite) );
01458
01459 std::multimap< unsigned, std::pair<double, unsigned> > associatedHOs;
01460
01461 PFClusterRef hclusterref = elements[iHcal].clusterRef();
01462 assert(!hclusterref.isNull() );
01463
01464
01465
01466
01467
01468
01469
01470
01471
01472 if( sortedTracks.empty() ) {
01473 if(debug_)
01474 cout<<"\tno associated tracks, keep for later"<<endl;
01475 continue;
01476 }
01477
01478
01479 active[iHcal] = false;
01480
01481
01482
01483
01484
01485
01486
01487 if(debug_) cout<<"\t"<<sortedTracks.size()<<" associated tracks:"<<endl;
01488
01489 double totalChargedMomentum = 0;
01490 double sumpError2 = 0.;
01491 double totalHO = 0.;
01492 double totalEcal = 0.;
01493 double totalHcal = hclusterref->energy();
01494 vector<double> hcalP;
01495 vector<double> hcalDP;
01496 vector<unsigned> tkIs;
01497 double maxDPovP = -9999.;
01498
01499
01500 vector< unsigned > chargedHadronsIndices;
01501 vector< unsigned > chargedHadronsInBlock;
01502 double mergedNeutralHadronEnergy = 0;
01503 double mergedPhotonEnergy = 0;
01504 double muonHCALEnergy = 0.;
01505 double muonECALEnergy = 0.;
01506 double muonHCALError = 0.;
01507 double muonECALError = 0.;
01508 unsigned nMuons = 0;
01509
01510
01511 math::XYZVector photonAtECAL(0.,0.,0.);
01512 math::XYZVector hadronDirection(hclusterref->position().X(),
01513 hclusterref->position().Y(),
01514 hclusterref->position().Z());
01515 hadronDirection = hadronDirection.Unit();
01516 math::XYZVector hadronAtECAL = totalHcal * hadronDirection;
01517
01518
01519 for(IE ie = sortedTracks.begin(); ie != sortedTracks.end(); ++ie ) {
01520
01521 unsigned iTrack = ie->second;
01522
01523
01524 if ( !active[iTrack] ) continue;
01525
01526 PFBlockElement::Type type = elements[iTrack].type();
01527 assert( type == reco::PFBlockElement::TRACK );
01528
01529 reco::TrackRef trackRef = elements[iTrack].trackRef();
01530 assert( !trackRef.isNull() );
01531
01532
01533 const math::XYZPointF& chargedPosition =
01534 dynamic_cast<const reco::PFBlockElementTrack*>(&elements[iTrack])->positionAtECALEntrance();
01535 math::XYZVector chargedDirection(chargedPosition.X(),chargedPosition.Y(),chargedPosition.Z());
01536 chargedDirection = chargedDirection.Unit();
01537
01538
01539 std::multimap<double, unsigned> sortedEcals;
01540 block.associatedElements( iTrack, linkData,
01541 sortedEcals,
01542 reco::PFBlockElement::ECAL,
01543 reco::PFBlock::LINKTEST_ALL );
01544
01545 if(debug_) cout<<"\t\t\tnumber of Ecal elements linked to this track: "
01546 <<sortedEcals.size()<<endl;
01547
01548
01549 std::multimap<double, unsigned> sortedHOs;
01550 if (useHO_) {
01551 block.associatedElements( iTrack, linkData,
01552 sortedHOs,
01553 reco::PFBlockElement::HO,
01554 reco::PFBlock::LINKTEST_ALL );
01555
01556 }
01557 if(debug_)
01558 cout<<"PFAlgo : number of HO elements linked to this track: "
01559 <<sortedHOs.size()<<endl;
01560
01561
01562 reco::MuonRef muonRef = elements[iTrack].muonRef();
01563
01564
01565 bool thisIsAMuon = PFMuonAlgo::isMuon(elements[iTrack]);
01566 bool thisIsAnIsolatedMuon = PFMuonAlgo::isIsolatedMuon(elements[iTrack]);
01567 bool thisIsALooseMuon = false;
01568
01569
01570 if(!thisIsAMuon ) {
01571 thisIsALooseMuon = PFMuonAlgo::isLooseMuon(elements[iTrack]);
01572 }
01573
01574
01575 if ( thisIsAMuon ) {
01576 if ( debug_ ) {
01577 std::cout << "\t\tThis track is identified as a muon - remove it from the stack" << std::endl;
01578 std::cout << "\t\t" << elements[iTrack] << std::endl;
01579 }
01580
01581
01582
01583 unsigned tmpi = reconstructTrack( elements[iTrack] );
01584
01585
01586 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iTrack );
01587 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iHcal );
01588 double muonHcal = std::min(muonHCAL_[0]+muonHCAL_[1],totalHcal);
01589
01590
01591
01592
01593 bool letMuonEatCaloEnergy = false;
01594
01595 if(thisIsAnIsolatedMuon){
01596
01597 double totalCaloEnergy = totalHcal / 1.30;
01598 unsigned iEcal = 0;
01599 if( !sortedEcals.empty() ) {
01600 iEcal = sortedEcals.begin()->second;
01601 PFClusterRef eclusterref = elements[iEcal].clusterRef();
01602 totalCaloEnergy += eclusterref->energy();
01603 }
01604
01605 if (useHO_) {
01606
01607 unsigned iHO = 0;
01608 if( !sortedHOs.empty() ) {
01609 iHO = sortedHOs.begin()->second;
01610 PFClusterRef eclusterref = elements[iHO].clusterRef();
01611 totalCaloEnergy += eclusterref->energy() / 1.30;
01612 }
01613 }
01614
01615
01616
01617 if( (pfCandidates_->back()).p() > totalCaloEnergy ) letMuonEatCaloEnergy = true;
01618 }
01619
01620 if(letMuonEatCaloEnergy) muonHcal = totalHcal;
01621 double muonEcal =0.;
01622 unsigned iEcal = 0;
01623 if( !sortedEcals.empty() ) {
01624 iEcal = sortedEcals.begin()->second;
01625 PFClusterRef eclusterref = elements[iEcal].clusterRef();
01626 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iEcal);
01627 muonEcal = std::min(muonECAL_[0]+muonECAL_[1],eclusterref->energy());
01628 if(letMuonEatCaloEnergy) muonEcal = eclusterref->energy();
01629
01630 if ( eclusterref->energy() - muonEcal < 0.2 ) active[iEcal] = false;
01631 (*pfCandidates_)[tmpi].setEcalEnergy(eclusterref->energy(), muonEcal);
01632 }
01633 unsigned iHO = 0;
01634 double muonHO =0.;
01635 if (useHO_) {
01636 if( !sortedHOs.empty() ) {
01637 iHO = sortedHOs.begin()->second;
01638 PFClusterRef hoclusterref = elements[iHO].clusterRef();
01639 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iHO);
01640 muonHO = std::min(muonHO_[0]+muonHO_[1],hoclusterref->energy());
01641 if(letMuonEatCaloEnergy) muonHO = hoclusterref->energy();
01642
01643 if ( hoclusterref->energy() - muonHO < 0.2 ) active[iHO] = false;
01644 (*pfCandidates_)[tmpi].setHcalEnergy(totalHcal, muonHcal);
01645 (*pfCandidates_)[tmpi].setHoEnergy(hoclusterref->energy(), muonHO);
01646 }
01647 } else {
01648 (*pfCandidates_)[tmpi].setHcalEnergy(totalHcal, muonHcal);
01649 }
01650
01651 if(letMuonEatCaloEnergy){
01652 muonHCALEnergy += totalHcal;
01653 if (useHO_) muonHCALEnergy +=muonHO;
01654 muonHCALError += 0.;
01655 muonECALEnergy += muonEcal;
01656 muonECALError += 0.;
01657 photonAtECAL -= muonEcal*chargedDirection;
01658 hadronAtECAL -= totalHcal*chargedDirection;
01659 if ( !sortedEcals.empty() ) active[iEcal] = false;
01660 active[iHcal] = false;
01661 if (useHO_ && !sortedHOs.empty() ) active[iHO] = false;
01662 }
01663 else{
01664
01665 muonHCALEnergy += muonHCAL_[0];
01666 muonHCALError += muonHCAL_[1]*muonHCAL_[1];
01667 if ( muonHO > 0. ) {
01668 muonHCALEnergy += muonHO_[0];
01669 muonHCALError += muonHO_[1]*muonHO_[1];
01670 }
01671 muonECALEnergy += muonECAL_[0];
01672 muonECALError += muonECAL_[1]*muonECAL_[1];
01673
01674 photonAtECAL -= muonECAL_[0]*chargedDirection;
01675 hadronAtECAL -= muonHCAL_[0]*chargedDirection;
01676 }
01677
01678
01679 active[iTrack] = false;
01680
01681 continue;
01682 }
01683
01684
01685
01686 if(debug_) cout<<"\t\t"<<elements[iTrack]<<endl;
01687
01688
01689 double trackMomentum = trackRef->p();
01690 totalChargedMomentum += trackMomentum;
01691
01692
01693
01694 if ( thisIsALooseMuon && !thisIsAMuon ) nMuons += 1;
01695
01696
01697
01698
01699 double Dpt = trackRef->ptError();
01700 double blowError = 1.;
01701 switch (trackRef->algo()) {
01702 case TrackBase::ctf:
01703 case TrackBase::iter0:
01704 case TrackBase::iter1:
01705 case TrackBase::iter2:
01706 case TrackBase::iter3:
01707 case TrackBase::iter4:
01708 blowError = 1.;
01709 break;
01710 case TrackBase::iter5:
01711 blowError = factors45_[0];
01712 break;
01713 case TrackBase::iter6:
01714 blowError = factors45_[1];
01715 break;
01716 default:
01717 blowError = 1E9;
01718 break;
01719 }
01720
01721 bool isPrimaryOrSecondary = isFromSecInt(elements[iTrack], "all");
01722
01723 if ( isPrimaryOrSecondary ) blowError = 1.;
01724
01725 std::pair<unsigned,bool> tkmuon = make_pair(iTrack,thisIsALooseMuon);
01726 associatedTracks.insert( make_pair(-Dpt*blowError, tkmuon) );
01727
01728
01729 double Dp = trackRef->qoverpError()*trackMomentum*trackMomentum;
01730 sumpError2 += Dp*Dp;
01731
01732 bool connectedToEcal = false;
01733 if( !sortedEcals.empty() )
01734 {
01735
01736
01737 for ( IE iec=sortedEcals.begin();
01738 iec!=sortedEcals.end(); ++iec ) {
01739
01740 unsigned iEcal = iec->second;
01741 double dist = iec->first;
01742
01743
01744 if( !active[iEcal] ) {
01745 if(debug_) cout<<"\t\tcluster locked"<<endl;
01746 continue;
01747 }
01748
01749
01750 PFBlockElement::Type type = elements[ iEcal ].type();
01751 assert( type == PFBlockElement::ECAL );
01752 PFClusterRef eclusterref = elements[iEcal].clusterRef();
01753 assert(!eclusterref.isNull() );
01754
01755
01756 std::multimap<double, unsigned> sortedTracksEcal;
01757 block.associatedElements( iEcal, linkData,
01758 sortedTracksEcal,
01759 reco::PFBlockElement::TRACK,
01760 reco::PFBlock::LINKTEST_ALL );
01761 unsigned jTrack = sortedTracksEcal.begin()->second;
01762 if ( jTrack != iTrack ) continue;
01763
01764
01765
01766 double distEcal = block.dist(jTrack,iEcal,linkData,
01767 reco::PFBlock::LINKTEST_ALL);
01768
01769
01770
01771
01772
01773
01774 vector<double> ps1Ene(1,static_cast<double>(0.));
01775 associatePSClusters(iEcal, reco::PFBlockElement::PS1,
01776 block, elements, linkData, active,
01777 ps1Ene);
01778 vector<double> ps2Ene(1,static_cast<double>(0.));
01779 associatePSClusters(iEcal, reco::PFBlockElement::PS2,
01780 block, elements, linkData, active,
01781 ps2Ene);
01782 std::pair<double,double> psEnes = make_pair(ps1Ene[0],
01783 ps2Ene[0]);
01784 associatedPSs.insert(make_pair(iEcal,psEnes));
01785
01786
01787 bool crackCorrection = false;
01788 float ecalEnergyCalibrated = calibration_->energyEm(*eclusterref,ps1Ene,ps2Ene,crackCorrection);
01789 math::XYZVector photonDirection(eclusterref->position().X(),
01790 eclusterref->position().Y(),
01791 eclusterref->position().Z());
01792 photonDirection = photonDirection.Unit();
01793
01794 if ( !connectedToEcal ) {
01795
01796 if(debug_) cout<<"\t\t\tclosest: "
01797 <<elements[iEcal]<<endl;
01798
01799 connectedToEcal = true;
01800
01801
01802
01803 std::pair<unsigned,math::XYZVector> satellite =
01804 make_pair(iEcal,ecalEnergyCalibrated*photonDirection);
01805 ecalSatellites.insert( make_pair(-1., satellite) );
01806
01807 } else {
01808
01809 std::pair<unsigned,math::XYZVector> satellite =
01810 make_pair(iEcal,ecalEnergyCalibrated*photonDirection);
01811 ecalSatellites.insert( make_pair(dist, satellite) );
01812
01813 }
01814
01815 std::pair<double, unsigned> associatedEcal
01816 = make_pair( distEcal, iEcal );
01817 associatedEcals.insert( make_pair(iTrack, associatedEcal) );
01818
01819 }
01820 }
01821
01822 if( useHO_ && !sortedHOs.empty() )
01823 {
01824
01825
01826 for ( IE ieho=sortedHOs.begin(); ieho!=sortedHOs.end(); ++ieho ) {
01827
01828 unsigned iHO = ieho->second;
01829 double distHO = ieho->first;
01830
01831
01832 if( !active[iHO] ) {
01833 if(debug_) cout<<"\t\tcluster locked"<<endl;
01834 continue;
01835 }
01836
01837
01838 PFBlockElement::Type type = elements[ iHO ].type();
01839 assert( type == PFBlockElement::HO );
01840 PFClusterRef hoclusterref = elements[iHO].clusterRef();
01841 assert(!hoclusterref.isNull() );
01842
01843
01844 std::multimap<double, unsigned> sortedTracksHO;
01845 block.associatedElements( iHO, linkData,
01846 sortedTracksHO,
01847 reco::PFBlockElement::TRACK,
01848 reco::PFBlock::LINKTEST_ALL );
01849 unsigned jTrack = sortedTracksHO.begin()->second;
01850 if ( jTrack != iTrack ) continue;
01851
01852
01853
01854
01855
01856
01857
01858 totalHO += hoclusterref->energy();
01859 active[iHO] = false;
01860
01861 std::pair<double, unsigned> associatedHO
01862 = make_pair( distHO, iHO );
01863 associatedHOs.insert( make_pair(iTrack, associatedHO) );
01864
01865 }
01866 }
01867
01868
01869 }
01870
01871
01872 totalHcal += totalHO;
01873
01874
01875
01876 double caloEnergy = 0.;
01877 double slopeEcal = 1.0;
01878 double calibEcal = 0.;
01879 double calibHcal = 0.;
01880 hadronDirection = hadronAtECAL.Unit();
01881
01882
01883 double Caloresolution = neutralHadronEnergyResolution( totalChargedMomentum, hclusterref->positionREP().Eta());
01884 Caloresolution *= totalChargedMomentum;
01885
01886 Caloresolution = std::sqrt(Caloresolution*Caloresolution + muonHCALError + muonECALError);
01887 totalEcal -= std::min(totalEcal,muonECALEnergy);
01888 totalHcal -= std::min(totalHcal,muonHCALEnergy);
01889 if ( totalEcal < 1E-9 ) photonAtECAL = math::XYZVector(0.,0.,0.);
01890 if ( totalHcal < 1E-9 ) hadronAtECAL = math::XYZVector(0.,0.,0.);
01891
01892
01893
01894
01895
01896 for ( IS is = ecalSatellites.begin(); is != ecalSatellites.end(); ++is ) {
01897
01898
01899 double previousCalibEcal = calibEcal;
01900 double previousCalibHcal = calibHcal;
01901 double previousCaloEnergy = caloEnergy;
01902 double previousSlopeEcal = slopeEcal;
01903 math::XYZVector previousHadronAtECAL = hadronAtECAL;
01904
01905 totalEcal += sqrt(is->second.second.Mag2());
01906 photonAtECAL += is->second.second;
01907 calibEcal = std::max(0.,totalEcal);
01908 calibHcal = std::max(0.,totalHcal);
01909 hadronAtECAL = calibHcal * hadronDirection;
01910
01911 calibration_->energyEmHad(totalChargedMomentum,calibEcal,calibHcal,
01912 hclusterref->positionREP().Eta(),
01913 hclusterref->positionREP().Phi());
01914 caloEnergy = calibEcal+calibHcal;
01915 if ( totalEcal > 0.) slopeEcal = calibEcal/totalEcal;
01916
01917 hadronAtECAL = calibHcal * hadronDirection;
01918
01919
01920
01921 if ( is->first < 0. || caloEnergy - totalChargedMomentum <= 0. ) {
01922 if(debug_) cout<<"\t\t\tactive, adding "<<is->second.second
01923 <<" to ECAL energy, and locking"<<endl;
01924 active[is->second.first] = false;
01925 continue;
01926 }
01927
01928
01929
01930 totalEcal -= sqrt(is->second.second.Mag2());
01931 photonAtECAL -= is->second.second;
01932 calibEcal = previousCalibEcal;
01933 calibHcal = previousCalibHcal;
01934 hadronAtECAL = previousHadronAtECAL;
01935 slopeEcal = previousSlopeEcal;
01936 caloEnergy = previousCaloEnergy;
01937
01938 break;
01939
01940 }
01941
01942
01943 assert(caloEnergy>=0);
01944
01945
01946
01947
01948
01949
01950
01951
01952
01953
01954
01955
01956 double TotalError = sqrt(sumpError2 + Caloresolution*Caloresolution);
01957
01958
01960 if ( totalChargedMomentum - caloEnergy > nSigmaTRACK_*Caloresolution ) {
01961
01962
01963
01964
01965
01966
01967
01968
01969
01970
01971
01972
01973
01974 if ( nMuons > 0 ) {
01975 for ( IT it = associatedTracks.begin(); it != associatedTracks.end(); ++it ) {
01976 unsigned iTrack = it->second.first;
01977
01978 if ( !active[iTrack] ) continue;
01979
01980 if ( !it->second.second ) continue;
01981
01982 double trackMomentum = elements[it->second.first].trackRef()->p();
01983
01984
01985 std::multimap<double, unsigned> sortedEcals;
01986 block.associatedElements( iTrack, linkData,
01987 sortedEcals,
01988 reco::PFBlockElement::ECAL,
01989 reco::PFBlock::LINKTEST_ALL );
01990 std::multimap<double, unsigned> sortedHOs;
01991 block.associatedElements( iTrack, linkData,
01992 sortedHOs,
01993 reco::PFBlockElement::HO,
01994 reco::PFBlock::LINKTEST_ALL );
01995
01996
01997 unsigned tmpi = reconstructTrack( elements[iTrack],true);
01998
01999 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iTrack );
02000 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iHcal );
02001 double muonHcal = std::min(muonHCAL_[0]+muonHCAL_[1],totalHcal-totalHO);
02002 double muonHO = 0.;
02003 (*pfCandidates_)[tmpi].setHcalEnergy(totalHcal,muonHcal);
02004 if( !sortedEcals.empty() ) {
02005 unsigned iEcal = sortedEcals.begin()->second;
02006 PFClusterRef eclusterref = elements[iEcal].clusterRef();
02007 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iEcal);
02008 double muonEcal = std::min(muonECAL_[0]+muonECAL_[1],eclusterref->energy());
02009 (*pfCandidates_)[tmpi].setEcalEnergy(eclusterref->energy(),muonEcal);
02010 }
02011 if( useHO_ && !sortedHOs.empty() ) {
02012 unsigned iHO = sortedHOs.begin()->second;
02013 PFClusterRef hoclusterref = elements[iHO].clusterRef();
02014 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iHO);
02015 muonHO = std::min(muonHO_[0]+muonHO_[1],hoclusterref->energy());
02016 (*pfCandidates_)[tmpi].setHcalEnergy(totalHcal-totalHO,muonHcal);
02017 (*pfCandidates_)[tmpi].setHoEnergy(hoclusterref->energy(),muonHO);
02018 }
02019
02020 const math::XYZPointF& chargedPosition =
02021 dynamic_cast<const reco::PFBlockElementTrack*>(&elements[it->second.first])->positionAtECALEntrance();
02022 math::XYZVector chargedDirection(chargedPosition.X(), chargedPosition.Y(), chargedPosition.Z());
02023 chargedDirection = chargedDirection.Unit();
02024 totalChargedMomentum -= trackMomentum;
02025
02026 if ( totalEcal > 0. )
02027 calibEcal -= std::min(calibEcal,muonECAL_[0]*calibEcal/totalEcal);
02028 if ( totalHcal > 0. )
02029 calibHcal -= std::min(calibHcal,muonHCAL_[0]*calibHcal/totalHcal);
02030 totalEcal -= std::min(totalEcal,muonECAL_[0]);
02031 totalHcal -= std::min(totalHcal,muonHCAL_[0]);
02032 if ( totalEcal > muonECAL_[0] ) photonAtECAL -= muonECAL_[0] * chargedDirection;
02033 if ( totalHcal > muonHCAL_[0] ) hadronAtECAL -= muonHCAL_[0]*calibHcal/totalHcal * chargedDirection;
02034 caloEnergy = calibEcal+calibHcal;
02035 muonHCALEnergy += muonHCAL_[0];
02036 muonHCALError += muonHCAL_[1]*muonHCAL_[1];
02037 if ( muonHO > 0. ) {
02038 muonHCALEnergy += muonHO_[0];
02039 muonHCALError += muonHO_[1]*muonHO_[1];
02040 if ( totalHcal > 0. ) {
02041 calibHcal -= std::min(calibHcal,muonHO_[0]*calibHcal/totalHcal);
02042 totalHcal -= std::min(totalHcal,muonHO_[0]);
02043 }
02044 }
02045 muonECALEnergy += muonECAL_[0];
02046 muonECALError += muonECAL_[1]*muonECAL_[1];
02047 active[iTrack] = false;
02048
02049
02050
02051 }
02052
02053 Caloresolution = neutralHadronEnergyResolution( totalChargedMomentum, hclusterref->positionREP().Eta());
02054 Caloresolution *= totalChargedMomentum;
02055 Caloresolution = std::sqrt(Caloresolution*Caloresolution + muonHCALError + muonECALError);
02056 }
02057 }
02058
02059
02060
02061
02062
02063
02064
02065
02066
02067
02068
02069
02070
02071
02072
02073 unsigned corrTrack = 10000000;
02074 double corrFact = 1.;
02075
02076 if (rejectTracks_Bad_ &&
02077 totalChargedMomentum - caloEnergy > nSigmaTRACK_*Caloresolution) {
02078
02079 for ( IT it = associatedTracks.begin(); it != associatedTracks.end(); ++it ) {
02080 unsigned iTrack = it->second.first;
02081
02082 if ( !active[iTrack] ) continue;
02083 const reco::TrackRef& trackref = elements[it->second.first].trackRef();
02084
02085 double dptRel = fabs(it->first)/trackref->pt()*100;
02086 bool isSecondary = isFromSecInt(elements[iTrack], "secondary");
02087 bool isPrimary = isFromSecInt(elements[iTrack], "primary");
02088
02089 if ( isSecondary && dptRel < dptRel_DispVtx_ ) continue;
02090
02091 if ( fabs(it->first) < ptError_ ) break;
02092
02093 double wouldBeTotalChargedMomentum =
02094 totalChargedMomentum - trackref->p();
02095
02096
02097
02098 if ( wouldBeTotalChargedMomentum > caloEnergy ) {
02099
02100 if (debug_ && isSecondary) {
02101 cout << "In bad track rejection step dptRel = " << dptRel << " dptRel_DispVtx_ = " << dptRel_DispVtx_ << endl;
02102 cout << "The calo energy would be still smaller even without this track but it is attached to a NI"<< endl;
02103 }
02104
02105
02106 if(isPrimary || (isSecondary && dptRel < dptRel_DispVtx_)) continue;
02107 active[iTrack] = false;
02108 totalChargedMomentum = wouldBeTotalChargedMomentum;
02109 if ( debug_ )
02110 std::cout << "\tElement " << elements[iTrack]
02111 << " rejected (Dpt = " << -it->first
02112 << " GeV/c, algo = " << trackref->algo() << ")" << std::endl;
02113
02114 } else {
02115 if(isPrimary) break;
02116 corrTrack = iTrack;
02117 corrFact = (caloEnergy - wouldBeTotalChargedMomentum)/elements[it->second.first].trackRef()->p();
02118 if ( trackref->p()*corrFact < 0.05 ) {
02119 corrFact = 0.;
02120 active[iTrack] = false;
02121 }
02122 totalChargedMomentum -= trackref->p()*(1.-corrFact);
02123 if ( debug_ )
02124 std::cout << "\tElement " << elements[iTrack]
02125 << " (Dpt = " << -it->first
02126 << " GeV/c, algo = " << trackref->algo()
02127 << ") rescaled by " << corrFact
02128 << " Now the total charged momentum is " << totalChargedMomentum << endl;
02129 break;
02130 }
02131 }
02132 }
02133
02134
02135 Caloresolution = neutralHadronEnergyResolution( totalChargedMomentum, hclusterref->positionREP().Eta());
02136 Caloresolution *= totalChargedMomentum;
02137 Caloresolution = std::sqrt(Caloresolution*Caloresolution + muonHCALError + muonECALError);
02138
02139
02140
02141
02142 if ( rejectTracks_Step45_ &&
02143 sortedTracks.size() > 1 &&
02144 totalChargedMomentum - caloEnergy > nSigmaTRACK_*Caloresolution ) {
02145 for ( IT it = associatedTracks.begin(); it != associatedTracks.end(); ++it ) {
02146 unsigned iTrack = it->second.first;
02147 reco::TrackRef trackref = elements[iTrack].trackRef();
02148 if ( !active[iTrack] ) continue;
02149
02150 double dptRel = fabs(it->first)/trackref->pt()*100;
02151 bool isPrimaryOrSecondary = isFromSecInt(elements[iTrack], "all");
02152
02153
02154
02155
02156 if (isPrimaryOrSecondary && dptRel < dptRel_DispVtx_) continue;
02157
02158 switch (trackref->algo()) {
02159 case TrackBase::ctf:
02160 case TrackBase::iter0:
02161 case TrackBase::iter1:
02162 case TrackBase::iter2:
02163 case TrackBase::iter3:
02164 case TrackBase::iter4:
02165 break;
02166 case TrackBase::iter5:
02167 case TrackBase::iter6:
02168 active[iTrack] = false;
02169 totalChargedMomentum -= trackref->p();
02170
02171 if ( debug_ )
02172 std::cout << "\tElement " << elements[iTrack]
02173 << " rejected (Dpt = " << -it->first
02174 << " GeV/c, algo = " << trackref->algo() << ")" << std::endl;
02175 break;
02176 default:
02177 break;
02178 }
02179 }
02180 }
02181
02182
02183 Caloresolution = neutralHadronEnergyResolution( totalChargedMomentum, hclusterref->positionREP().Eta());
02184 Caloresolution *= totalChargedMomentum;
02185 Caloresolution = std::sqrt(Caloresolution*Caloresolution + muonHCALError + muonECALError);
02186
02187
02188 sumpError2 = 0.;
02189 for ( IT it = associatedTracks.begin(); it != associatedTracks.end(); ++it ) {
02190 unsigned iTrack = it->second.first;
02191 if ( !active[iTrack] ) continue;
02192 reco::TrackRef trackRef = elements[iTrack].trackRef();
02193 double trackMomentum = trackRef->p();
02194 double Dp = trackRef->qoverpError()*trackMomentum*trackMomentum;
02195 unsigned tmpi = reconstructTrack( elements[iTrack] );
02196
02197
02198 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iTrack );
02199 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iHcal );
02200 std::pair<II,II> myEcals = associatedEcals.equal_range(iTrack);
02201 for (II ii=myEcals.first; ii!=myEcals.second; ++ii ) {
02202 unsigned iEcal = ii->second.second;
02203 if ( active[iEcal] ) continue;
02204 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iEcal );
02205 }
02206
02207 if (useHO_) {
02208 std::pair<II,II> myHOs = associatedHOs.equal_range(iTrack);
02209 for (II ii=myHOs.first; ii!=myHOs.second; ++ii ) {
02210 unsigned iHO = ii->second.second;
02211 if ( active[iHO] ) continue;
02212 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iHO );
02213 }
02214 }
02215
02216 if ( iTrack == corrTrack ) {
02217 (*pfCandidates_)[tmpi].rescaleMomentum(corrFact);
02218 trackMomentum *= corrFact;
02219 }
02220 chargedHadronsIndices.push_back( tmpi );
02221 chargedHadronsInBlock.push_back( iTrack );
02222 active[iTrack] = false;
02223 hcalP.push_back(trackMomentum);
02224 hcalDP.push_back(Dp);
02225 if (Dp/trackMomentum > maxDPovP) maxDPovP = Dp/trackMomentum;
02226 sumpError2 += Dp*Dp;
02227 }
02228
02229
02230 TotalError = sqrt(sumpError2 + Caloresolution*Caloresolution);
02231
02232 if ( debug_ ) {
02233 cout<<"\tCompare Calo Energy to total charged momentum "<<endl;
02234 cout<<"\t\tsum p = "<<totalChargedMomentum<<" +- "<<sqrt(sumpError2)<<endl;
02235 cout<<"\t\tsum ecal = "<<totalEcal<<endl;
02236 cout<<"\t\tsum hcal = "<<totalHcal<<endl;
02237 cout<<"\t\t => Calo Energy = "<<caloEnergy<<" +- "<<Caloresolution<<endl;
02238 cout<<"\t\t => Calo Energy- total charged momentum = "
02239 <<caloEnergy-totalChargedMomentum<<" +- "<<TotalError<<endl;
02240 cout<<endl;
02241 }
02242
02243
02244
02246 double nsigma = nSigmaHCAL(totalChargedMomentum,hclusterref->positionREP().Eta());
02247
02248 if ( abs(totalChargedMomentum-caloEnergy)<nsigma*TotalError ) {
02249
02250
02251
02252
02253 if(debug_) {
02254 cout<<"\t\tcase 1: COMPATIBLE "
02255 <<"|Calo Energy- total charged momentum| = "
02256 <<abs(caloEnergy-totalChargedMomentum)
02257 <<" < "<<nsigma<<" x "<<TotalError<<endl;
02258 if (maxDPovP < 0.1 )
02259 cout<<"\t\t\tmax DP/P = "<< maxDPovP
02260 <<" less than 0.1: do nothing "<<endl;
02261 else
02262 cout<<"\t\t\tmax DP/P = "<< maxDPovP
02263 <<" > 0.1: take weighted averages "<<endl;
02264 }
02265
02266
02267 if (maxDPovP > 0.1) {
02268
02269
02270
02271 int nrows = chargedHadronsIndices.size();
02272 TMatrixTSym<double> a (nrows);
02273 TVectorD b (nrows);
02274 TVectorD check(nrows);
02275 double sigma2E = Caloresolution*Caloresolution;
02276 for(int i=0; i<nrows; i++) {
02277 double sigma2i = hcalDP[i]*hcalDP[i];
02278 if (debug_){
02279 cout<<"\t\t\ttrack associated to hcal "<<i
02280 <<" P = "<<hcalP[i]<<" +- "
02281 <<hcalDP[i]<<endl;
02282 }
02283 a(i,i) = 1./sigma2i + 1./sigma2E;
02284 b(i) = hcalP[i]/sigma2i + caloEnergy/sigma2E;
02285 for(int j=0; j<nrows; j++) {
02286 if (i==j) continue;
02287 a(i,j) = 1./sigma2E;
02288 }
02289 }
02290
02291
02292
02293
02294
02295
02296
02297
02298 TDecompChol decomp(a);
02299 bool ok = false;
02300 TVectorD x = decomp.Solve(b, ok);
02301
02302
02303 if (ok) {
02304 for (int i=0; i<nrows; i++){
02305
02306 unsigned ich = chargedHadronsIndices[i];
02307 double rescaleFactor = x(i)/hcalP[i];
02308 (*pfCandidates_)[ich].rescaleMomentum( rescaleFactor );
02309
02310 if(debug_){
02311 cout<<"\t\t\told p "<<hcalP[i]
02312 <<" new p "<<x(i)
02313 <<" rescale "<<rescaleFactor<<endl;
02314 }
02315 }
02316 }
02317 else {
02318 cerr<<"not ok!"<<endl;
02319 assert(0);
02320 }
02321 }
02322 }
02323
02325 else if( caloEnergy > totalChargedMomentum ) {
02326
02327
02328
02329
02330
02331
02332
02333
02334
02335
02336
02337
02338
02339
02340
02341
02342
02343
02344 double eNeutralHadron = caloEnergy - totalChargedMomentum;
02345 double ePhoton = (caloEnergy - totalChargedMomentum) / slopeEcal;
02346
02347 if(debug_) {
02348 if(!sortedTracks.empty() ){
02349 cout<<"\t\tcase 2: NEUTRAL DETECTION "
02350 <<caloEnergy<<" > "<<nsigma<<"x"<<TotalError
02351 <<" + "<<totalChargedMomentum<<endl;
02352 cout<<"\t\tneutral activity detected: "<<endl
02353 <<"\t\t\t photon = "<<ePhoton<<endl
02354 <<"\t\t\tor neutral hadron = "<<eNeutralHadron<<endl;
02355
02356 cout<<"\t\tphoton or hadron ?"<<endl;}
02357
02358 if(sortedTracks.empty() )
02359 cout<<"\t\tno track -> hadron "<<endl;
02360 else
02361 cout<<"\t\t"<<sortedTracks.size()
02362 <<"tracks -> check if the excess is photonic or hadronic"<<endl;
02363 }
02364
02365
02366 double ratioMax = 0.;
02367 reco::PFClusterRef maxEcalRef;
02368 unsigned maxiEcal= 9999;
02369
02370
02371
02372 for(IE ie = sortedTracks.begin(); ie != sortedTracks.end(); ++ie ) {
02373
02374 unsigned iTrack = ie->second;
02375
02376 PFBlockElement::Type type = elements[iTrack].type();
02377 assert( type == reco::PFBlockElement::TRACK );
02378
02379 reco::TrackRef trackRef = elements[iTrack].trackRef();
02380 assert( !trackRef.isNull() );
02381
02382 II iae = associatedEcals.find(iTrack);
02383
02384 if( iae == associatedEcals.end() ) continue;
02385
02386
02387 unsigned iEcal = iae->second.second;
02388
02389 PFBlockElement::Type typeEcal = elements[iEcal].type();
02390 assert(typeEcal == reco::PFBlockElement::ECAL);
02391
02392 reco::PFClusterRef clusterRef = elements[iEcal].clusterRef();
02393 assert( !clusterRef.isNull() );
02394
02395 double pTrack = trackRef->p();
02396 double eECAL = clusterRef->energy();
02397 double eECALOverpTrack = eECAL / pTrack;
02398
02399 if ( eECALOverpTrack > ratioMax ) {
02400 ratioMax = eECALOverpTrack;
02401 maxEcalRef = clusterRef;
02402 maxiEcal = iEcal;
02403 }
02404
02405 }
02406
02407 std::vector<reco::PFClusterRef> pivotalClusterRef;
02408 std::vector<unsigned> iPivotal;
02409 std::vector<double> particleEnergy, ecalEnergy, hcalEnergy, rawecalEnergy, rawhcalEnergy;
02410 std::vector<math::XYZVector> particleDirection;
02411
02412
02413
02414 if ( ePhoton < totalEcal || eNeutralHadron-calibEcal < 1E-10 ) {
02415
02416 if ( !maxEcalRef.isNull() ) {
02417 particleEnergy.push_back(ePhoton);
02418 particleDirection.push_back(photonAtECAL);
02419 ecalEnergy.push_back(ePhoton);
02420 hcalEnergy.push_back(0.);
02421 rawecalEnergy.push_back(totalEcal);
02422 rawhcalEnergy.push_back(totalHcal);
02423 pivotalClusterRef.push_back(maxEcalRef);
02424 iPivotal.push_back(maxiEcal);
02425
02426 mergedPhotonEnergy = ePhoton;
02427 }
02428
02429 } else {
02430
02431
02432 if ( !maxEcalRef.isNull() ) {
02433 pivotalClusterRef.push_back(maxEcalRef);
02434 particleEnergy.push_back(totalEcal);
02435 particleDirection.push_back(photonAtECAL);
02436 ecalEnergy.push_back(totalEcal);
02437 hcalEnergy.push_back(0.);
02438 rawecalEnergy.push_back(totalEcal);
02439 rawhcalEnergy.push_back(totalHcal);
02440 iPivotal.push_back(maxiEcal);
02441
02442 mergedPhotonEnergy = totalEcal;
02443 }
02444
02445
02446 mergedNeutralHadronEnergy = eNeutralHadron-calibEcal;
02447 if ( mergedNeutralHadronEnergy > 1.0 ) {
02448 pivotalClusterRef.push_back(hclusterref);
02449 particleEnergy.push_back(mergedNeutralHadronEnergy);
02450 particleDirection.push_back(hadronAtECAL);
02451 ecalEnergy.push_back(0.);
02452 hcalEnergy.push_back(mergedNeutralHadronEnergy);
02453 rawecalEnergy.push_back(totalEcal);
02454 rawhcalEnergy.push_back(totalHcal);
02455 iPivotal.push_back(iHcal);
02456 }
02457
02458 }
02459
02460
02461
02462
02463
02464
02465 for ( unsigned iPivot=0; iPivot<iPivotal.size(); ++iPivot ) {
02466
02467 if ( particleEnergy[iPivot] < 0. )
02468 std::cout << "ALARM = Negative energy ! "
02469 << particleEnergy[iPivot] << std::endl;
02470
02471 bool useDirection = true;
02472 unsigned tmpi = reconstructCluster( *pivotalClusterRef[iPivot],
02473 particleEnergy[iPivot],
02474 useDirection,
02475 particleDirection[iPivot].X(),
02476 particleDirection[iPivot].Y(),
02477 particleDirection[iPivot].Z());
02478
02479
02480 (*pfCandidates_)[tmpi].setEcalEnergy( rawecalEnergy[iPivot],ecalEnergy[iPivot] );
02481 if ( !useHO_ ) {
02482 (*pfCandidates_)[tmpi].setHcalEnergy( rawhcalEnergy[iPivot],hcalEnergy[iPivot] );
02483 (*pfCandidates_)[tmpi].setHoEnergy(0., 0.);
02484 } else {
02485 (*pfCandidates_)[tmpi].setHcalEnergy( rawhcalEnergy[iPivot]-totalHO,hcalEnergy[iPivot]*(1.-totalHO/rawhcalEnergy[iPivot]));
02486 (*pfCandidates_)[tmpi].setHoEnergy(totalHO, totalHO * hcalEnergy[iPivot]/rawhcalEnergy[iPivot]);
02487 }
02488 (*pfCandidates_)[tmpi].setPs1Energy( 0. );
02489 (*pfCandidates_)[tmpi].setPs2Energy( 0. );
02490 (*pfCandidates_)[tmpi].set_mva_nothing_gamma( -1. );
02491
02492
02493 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iHcal );
02494 for ( unsigned ich=0; ich<chargedHadronsInBlock.size(); ++ich) {
02495 unsigned iTrack = chargedHadronsInBlock[ich];
02496 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iTrack );
02497
02498 const math::XYZPointF& chargedPosition =
02499 dynamic_cast<const reco::PFBlockElementTrack*>(&elements[iTrack])->positionAtECALEntrance();
02500 (*pfCandidates_)[tmpi].setPositionAtECALEntrance(chargedPosition);
02501
02502 std::pair<II,II> myEcals = associatedEcals.equal_range(iTrack);
02503 for (II ii=myEcals.first; ii!=myEcals.second; ++ii ) {
02504 unsigned iEcal = ii->second.second;
02505 if ( active[iEcal] ) continue;
02506 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iEcal );
02507 }
02508 }
02509
02510 }
02511
02512 }
02513
02514
02515
02516
02517
02518 double totalHcalEnergyCalibrated = calibHcal;
02519 double totalEcalEnergyCalibrated = calibEcal;
02520
02521
02522
02523
02524
02525
02526
02527
02528
02529
02530
02531
02532
02533
02534 totalHcalEnergyCalibrated -= mergedNeutralHadronEnergy;
02535
02536 totalEcalEnergyCalibrated -= mergedPhotonEnergy;
02537
02538
02539
02540
02541 double chargedHadronsTotalEnergy = 0;
02542 for( unsigned ich=0; ich<chargedHadronsIndices.size(); ++ich ) {
02543 unsigned index = chargedHadronsIndices[ich];
02544 reco::PFCandidate& chargedHadron = (*pfCandidates_)[index];
02545 chargedHadronsTotalEnergy += chargedHadron.energy();
02546 }
02547
02548 for( unsigned ich=0; ich<chargedHadronsIndices.size(); ++ich ) {
02549 unsigned index = chargedHadronsIndices[ich];
02550 reco::PFCandidate& chargedHadron = (*pfCandidates_)[index];
02551 float fraction = chargedHadron.energy()/chargedHadronsTotalEnergy;
02552
02553 if ( !useHO_ ) {
02554 chargedHadron.setHcalEnergy( fraction * totalHcal, fraction * totalHcalEnergyCalibrated );
02555 chargedHadron.setHoEnergy( 0., 0. );
02556 } else {
02557 chargedHadron.setHcalEnergy( fraction * (totalHcal-totalHO), fraction * totalHcalEnergyCalibrated * (1.-totalHO/totalHcal) );
02558 chargedHadron.setHoEnergy( fraction * totalHO, fraction * totalHO * totalHcalEnergyCalibrated / totalHcal );
02559 }
02560
02561 chargedHadron.setEcalEnergy( fraction * totalEcal, fraction * totalEcalEnergyCalibrated );
02562 }
02563
02564
02565 for ( IS is = ecalSatellites.begin(); is != ecalSatellites.end(); ++is ) {
02566
02567
02568 unsigned iEcal = is->second.first;
02569 if ( !active[iEcal] ) continue;
02570
02571
02572 PFBlockElement::Type type = elements[ iEcal ].type();
02573 assert( type == PFBlockElement::ECAL );
02574 PFClusterRef eclusterref = elements[iEcal].clusterRef();
02575 assert(!eclusterref.isNull() );
02576
02577
02578 active[iEcal] = false;
02579
02580
02581 std::multimap<double, unsigned> assTracks;
02582 block.associatedElements( iEcal, linkData,
02583 assTracks,
02584 reco::PFBlockElement::TRACK,
02585 reco::PFBlock::LINKTEST_ALL );
02586
02587
02588 unsigned tmpi = reconstructCluster( *eclusterref, sqrt(is->second.second.Mag2()) );
02589 (*pfCandidates_)[tmpi].setEcalEnergy( eclusterref->energy(),sqrt(is->second.second.Mag2()) );
02590 (*pfCandidates_)[tmpi].setHcalEnergy( 0., 0. );
02591 (*pfCandidates_)[tmpi].setHoEnergy( 0., 0. );
02592 (*pfCandidates_)[tmpi].setPs1Energy( associatedPSs[iEcal].first );
02593 (*pfCandidates_)[tmpi].setPs2Energy( associatedPSs[iEcal].second );
02594 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iEcal );
02595 (*pfCandidates_)[tmpi].addElementInBlock( blockref, sortedTracks.begin()->second) ;
02596 }
02597
02598
02599 }
02600
02601
02602
02603
02604
02605 if(debug_) {
02606 cout<<endl;
02607 if(debug_)
02608 cout<<endl
02609 <<"---- loop remaining hcal ------- "<<endl;
02610 }
02611
02612
02613
02614
02615 for(unsigned ihcluster=0; ihcluster<hcalIs.size(); ihcluster++) {
02616 unsigned iHcal = hcalIs[ihcluster];
02617
02618
02619 std::vector<unsigned> ecalRefs;
02620 std::vector<unsigned> hoRefs;
02621
02622 if(debug_)
02623 cout<<endl<<elements[iHcal]<<" ";
02624
02625
02626 if( !active[iHcal] ) {
02627 if(debug_)
02628 cout<<"not active"<<endl;
02629 continue;
02630 }
02631
02632
02633 std::multimap<double, unsigned> ecalElems;
02634 block.associatedElements( iHcal, linkData,
02635 ecalElems ,
02636 reco::PFBlockElement::ECAL,
02637 reco::PFBlock::LINKTEST_ALL );
02638
02639
02640 float totalEcal = 0.;
02641 float ecalMax = 0.;
02642 reco::PFClusterRef eClusterRef;
02643 for(IE ie = ecalElems.begin(); ie != ecalElems.end(); ++ie ) {
02644
02645 unsigned iEcal = ie->second;
02646 double dist = ie->first;
02647 PFBlockElement::Type type = elements[iEcal].type();
02648 assert( type == PFBlockElement::ECAL );
02649
02650
02651 if( !active[iEcal] ) continue;
02652
02653
02654
02655
02656
02657
02658
02659
02660
02661
02662
02663
02664
02665
02666
02667
02668
02669
02670
02671
02672 std::multimap<double, unsigned> hcalElems;
02673 block.associatedElements( iEcal, linkData,
02674 hcalElems ,
02675 reco::PFBlockElement::HCAL,
02676 reco::PFBlock::LINKTEST_ALL );
02677
02678 bool isClosest = true;
02679 for(IE ih = hcalElems.begin(); ih != hcalElems.end(); ++ih ) {
02680
02681 unsigned jHcal = ih->second;
02682 double distH = ih->first;
02683
02684 if ( !active[jHcal] ) continue;
02685
02686 if ( distH < dist ) {
02687 isClosest = false;
02688 break;
02689 }
02690
02691 }
02692
02693 if (!isClosest) continue;
02694
02695
02696 if(debug_) {
02697 cout<<"\telement "<<elements[iEcal]<<" linked with dist "<< dist<<endl;
02698 cout<<"Added to HCAL cluster to form a neutral hadron"<<endl;
02699 }
02700
02701 reco::PFClusterRef eclusterRef = elements[iEcal].clusterRef();
02702 assert( !eclusterRef.isNull() );
02703
02704
02705 vector<double> ps1Ene(1,static_cast<double>(0.));
02706 associatePSClusters(iEcal, reco::PFBlockElement::PS1, block, elements, linkData, active, ps1Ene);
02707 vector<double> ps2Ene(1,static_cast<double>(0.));
02708 associatePSClusters(iEcal, reco::PFBlockElement::PS2, block, elements, linkData, active, ps2Ene);
02709 bool crackCorrection = false;
02710 double ecalEnergy = calibration_->energyEm(*eclusterRef,ps1Ene,ps2Ene,crackCorrection);
02711
02712
02713
02714 totalEcal += ecalEnergy;
02715 if ( ecalEnergy > ecalMax ) {
02716 ecalMax = ecalEnergy;
02717 eClusterRef = eclusterRef;
02718 }
02719
02720 ecalRefs.push_back(iEcal);
02721 active[iEcal] = false;
02722
02723
02724 }
02725
02726
02727 double totalHO = 0.;
02728 double hoMax = 0.;
02729
02730 if (useHO_) {
02731 std::multimap<double, unsigned> hoElems;
02732 block.associatedElements( iHcal, linkData,
02733 hoElems ,
02734 reco::PFBlockElement::HO,
02735 reco::PFBlock::LINKTEST_ALL );
02736
02737
02738
02739
02740
02741 reco::PFClusterRef hoClusterRef;
02742 for(IE ie = hoElems.begin(); ie != hoElems.end(); ++ie ) {
02743
02744 unsigned iHO = ie->second;
02745 double dist = ie->first;
02746 PFBlockElement::Type type = elements[iHO].type();
02747 assert( type == PFBlockElement::HO );
02748
02749
02750 if( !active[iHO] ) continue;
02751
02752
02753
02754
02755
02756 std::multimap<double, unsigned> hcalElems;
02757 block.associatedElements( iHO, linkData,
02758 hcalElems ,
02759 reco::PFBlockElement::HCAL,
02760 reco::PFBlock::LINKTEST_ALL );
02761
02762 bool isClosest = true;
02763 for(IE ih = hcalElems.begin(); ih != hcalElems.end(); ++ih ) {
02764
02765 unsigned jHcal = ih->second;
02766 double distH = ih->first;
02767
02768 if ( !active[jHcal] ) continue;
02769
02770 if ( distH < dist ) {
02771 isClosest = false;
02772 break;
02773 }
02774
02775 }
02776
02777 if (!isClosest) continue;
02778
02779 if(debug_ && useHO_) {
02780 cout<<"\telement "<<elements[iHO]<<" linked with dist "<< dist<<endl;
02781 cout<<"Added to HCAL cluster to form a neutral hadron"<<endl;
02782 }
02783
02784 reco::PFClusterRef hoclusterRef = elements[iHO].clusterRef();
02785 assert( !hoclusterRef.isNull() );
02786
02787 double hoEnergy = hoclusterRef->energy();
02788
02789 totalHO += hoEnergy;
02790 if ( hoEnergy > hoMax ) {
02791 hoMax = hoEnergy;
02792 hoClusterRef = hoclusterRef;
02793
02794 }
02795
02796 hoRefs.push_back(iHO);
02797 active[iHO] = false;
02798
02799 }
02800 }
02801
02802 PFClusterRef hclusterRef
02803 = elements[iHcal].clusterRef();
02804 assert( !hclusterRef.isNull() );
02805
02806
02807 double totalHcal = hclusterRef->energy();
02808
02809 if ( useHO_ ) totalHcal += totalHO;
02810
02811
02812
02813
02814
02815
02816
02817 double calibEcal = totalEcal > 0. ? totalEcal : 0.;
02818 double calibHcal = std::max(0.,totalHcal);
02819 if ( hclusterRef->layer() == PFLayer::HF_HAD ||
02820 hclusterRef->layer() == PFLayer::HF_EM ) {
02821
02822 calibEcal = totalEcal;
02823 } else {
02824 calibration_->energyEmHad(-1.,calibEcal,calibHcal,
02825 hclusterRef->positionREP().Eta(),
02826 hclusterRef->positionREP().Phi());
02827
02828 }
02829
02830
02831
02832
02833
02834
02835
02836
02837
02838 unsigned tmpi = reconstructCluster( *hclusterRef,
02839 calibEcal+calibHcal );
02840
02841
02842 (*pfCandidates_)[tmpi].setEcalEnergy( totalEcal, calibEcal );
02843 if ( !useHO_ ) {
02844 (*pfCandidates_)[tmpi].setHcalEnergy( totalHcal, calibHcal );
02845 (*pfCandidates_)[tmpi].setHoEnergy(0.,0.);
02846 } else {
02847 (*pfCandidates_)[tmpi].setHcalEnergy( totalHcal-totalHO, calibHcal*(1.-totalHO/totalHcal));
02848 (*pfCandidates_)[tmpi].setHoEnergy(totalHO,totalHO*calibHcal/totalHcal);
02849 }
02850 (*pfCandidates_)[tmpi].setPs1Energy( 0. );
02851 (*pfCandidates_)[tmpi].setPs2Energy( 0. );
02852 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iHcal );
02853 for (unsigned iec=0; iec<ecalRefs.size(); ++iec)
02854 (*pfCandidates_)[tmpi].addElementInBlock( blockref, ecalRefs[iec] );
02855 for (unsigned iho=0; iho<hoRefs.size(); ++iho)
02856 (*pfCandidates_)[tmpi].addElementInBlock( blockref, hoRefs[iho] );
02857
02858 }
02859
02860
02861
02862
02863 if(debug_) {
02864 cout<<endl;
02865 if(debug_) cout<<endl<<"---- loop ecal------- "<<endl;
02866 }
02867
02868
02869
02870 for(unsigned i=0; i<ecalIs.size(); i++) {
02871 unsigned iEcal = ecalIs[i];
02872
02873 if(debug_)
02874 cout<<endl<<elements[iEcal]<<" ";
02875
02876 if( ! active[iEcal] ) {
02877 if(debug_)
02878 cout<<"not active"<<endl;
02879 continue;
02880 }
02881
02882 PFBlockElement::Type type = elements[ iEcal ].type();
02883 assert( type == PFBlockElement::ECAL );
02884
02885 PFClusterRef clusterref = elements[iEcal].clusterRef();
02886 assert(!clusterref.isNull() );
02887
02888 active[iEcal] = false;
02889
02890
02891 vector<double> ps1Ene(1,static_cast<double>(0.));
02892 associatePSClusters(iEcal, reco::PFBlockElement::PS1, block, elements, linkData, active, ps1Ene);
02893 vector<double> ps2Ene(1,static_cast<double>(0.));
02894 associatePSClusters(iEcal, reco::PFBlockElement::PS2, block, elements, linkData, active, ps2Ene);
02895 bool crackCorrection = false;
02896 float ecalEnergy = calibration_->energyEm(*clusterref,ps1Ene,ps2Ene,crackCorrection);
02897
02898 double particleEnergy = ecalEnergy;
02899
02900 unsigned tmpi = reconstructCluster( *clusterref,
02901 particleEnergy );
02902
02903 (*pfCandidates_)[tmpi].setEcalEnergy( clusterref->energy(),ecalEnergy );
02904 (*pfCandidates_)[tmpi].setHcalEnergy( 0., 0. );
02905 (*pfCandidates_)[tmpi].setHoEnergy( 0., 0. );
02906 (*pfCandidates_)[tmpi].setPs1Energy( 0. );
02907 (*pfCandidates_)[tmpi].setPs2Energy( 0. );
02908 (*pfCandidates_)[tmpi].addElementInBlock( blockref, iEcal );
02909
02910
02911 }
02912
02913 }
02914
02916 unsigned PFAlgo::reconstructTrack( const reco::PFBlockElement& elt, bool allowLoose) {
02917
02918 const reco::PFBlockElementTrack* eltTrack
02919 = dynamic_cast<const reco::PFBlockElementTrack*>(&elt);
02920
02921 reco::TrackRef trackRef = eltTrack->trackRef();
02922 const reco::Track& track = *trackRef;
02923 reco::MuonRef muonRef = eltTrack->muonRef();
02924 int charge = track.charge()>0 ? 1 : -1;
02925
02926
02927 double px = track.px();
02928 double py = track.py();
02929 double pz = track.pz();
02930 double energy = sqrt(track.p()*track.p() + 0.13957*0.13957);
02931
02932
02933 math::XYZTLorentzVector momentum(px,py,pz,energy);
02934 reco::PFCandidate::ParticleType particleType
02935 = reco::PFCandidate::h;
02936
02937
02938 pfCandidates_->push_back( PFCandidate( charge,
02939 momentum,
02940 particleType ) );
02941
02942 pfCandidates_->back().setVertexSource( PFCandidate::kTrkVertex );
02943 pfCandidates_->back().setTrackRef( trackRef );
02944 pfCandidates_->back().setPositionAtECALEntrance( eltTrack->positionAtECALEntrance());
02945 if( muonRef.isNonnull())
02946 pfCandidates_->back().setMuonRef( muonRef );
02947
02948
02949
02950
02951 bool isMuon=pfmu_->reconstructMuon(pfCandidates_->back(),muonRef,allowLoose);
02952 bool isFromDisp = isFromSecInt(elt, "secondary");
02953
02954
02955 if ((!isMuon) && isFromDisp) {
02956 double Dpt = trackRef->ptError();
02957 double dptRel = Dpt/trackRef->pt()*100;
02958
02959
02960
02961 if (dptRel < dptRel_DispVtx_){
02962 if (debug_)
02963 cout << "Not refitted px = " << px << " py = " << py << " pz = " << pz << " energy = " << energy << endl;
02964
02965 reco::PFDisplacedVertexRef vRef = eltTrack->displacedVertexRef(reco::PFBlockElement::T_FROM_DISP)->displacedVertexRef();
02966 reco::Track trackRefit = vRef->refittedTrack(trackRef);
02967
02968 math::XYZTLorentzVector momentum(trackRefit.px(),
02969 trackRefit.py(),
02970 trackRefit.pz(),
02971 sqrt(trackRefit.p()*trackRefit.p() + 0.13957*0.13957));
02972 if (debug_)
02973 cout << "Refitted px = " << px << " py = " << py << " pz = " << pz << " energy = " << energy << endl;
02974 }
02975 pfCandidates_->back().setFlag( reco::PFCandidate::T_FROM_DISP, true);
02976 pfCandidates_->back().setDisplacedVertexRef( eltTrack->displacedVertexRef(reco::PFBlockElement::T_FROM_DISP)->displacedVertexRef(), reco::PFCandidate::T_FROM_DISP);
02977 }
02978
02979
02980 if(isFromSecInt(elt, "primary") && !isMuon) {
02981 pfCandidates_->back().setFlag( reco::PFCandidate::T_TO_DISP, true);
02982 pfCandidates_->back().setDisplacedVertexRef( eltTrack->displacedVertexRef(reco::PFBlockElement::T_TO_DISP)->displacedVertexRef(), reco::PFCandidate::T_TO_DISP);
02983 }
02984
02985 return pfCandidates_->size()-1;
02986 }
02987
02988
02989 unsigned
02990 PFAlgo::reconstructCluster(const reco::PFCluster& cluster,
02991 double particleEnergy,
02992 bool useDirection,
02993 double particleX,
02994 double particleY,
02995 double particleZ) {
02996
02997 reco::PFCandidate::ParticleType particleType = reco::PFCandidate::X;
02998
02999
03000
03001
03002
03003 double factor = 1.;
03004 if ( useDirection ) {
03005 switch( cluster.layer() ) {
03006 case PFLayer::ECAL_BARREL:
03007 case PFLayer::HCAL_BARREL1:
03008 factor = std::sqrt(cluster.position().Perp2()/(particleX*particleX+particleY*particleY));
03009 break;
03010 case PFLayer::ECAL_ENDCAP:
03011 case PFLayer::HCAL_ENDCAP:
03012 case PFLayer::HF_HAD:
03013 case PFLayer::HF_EM:
03014 factor = cluster.position().Z()/particleZ;
03015 break;
03016 default:
03017 assert(0);
03018 }
03019 }
03020
03021 math::XYZPoint vertexPos;
03022 if(useVertices_)
03023 vertexPos = math::XYZPoint(primaryVertex_.x(),primaryVertex_.y(),primaryVertex_.z());
03024 else
03025 vertexPos = math::XYZPoint(0.0,0.0,0.0);
03026
03027
03028 math::XYZVector clusterPos( cluster.position().X()-vertexPos.X(),
03029 cluster.position().Y()-vertexPos.Y(),
03030 cluster.position().Z()-vertexPos.Z());
03031 math::XYZVector particleDirection ( particleX*factor-vertexPos.X(),
03032 particleY*factor-vertexPos.Y(),
03033 particleZ*factor-vertexPos.Z() );
03034
03035
03036
03037
03038 clusterPos = useDirection ? particleDirection.Unit() : clusterPos.Unit();
03039 clusterPos *= particleEnergy;
03040
03041
03042
03043
03044 double mass = 0;
03045 ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> >
03046 momentum( clusterPos.X(), clusterPos.Y(), clusterPos.Z(), mass);
03047
03048 math::XYZTLorentzVector tmp;
03049
03050 tmp = momentum;
03051
03052
03053 int charge = 0;
03054
03055
03056 switch( cluster.layer() ) {
03057 case PFLayer::ECAL_BARREL:
03058 case PFLayer::ECAL_ENDCAP:
03059 particleType = PFCandidate::gamma;
03060 break;
03061 case PFLayer::HCAL_BARREL1:
03062 case PFLayer::HCAL_ENDCAP:
03063 particleType = PFCandidate::h0;
03064 break;
03065 case PFLayer::HF_HAD:
03066 particleType = PFCandidate::h_HF;
03067 break;
03068 case PFLayer::HF_EM:
03069 particleType = PFCandidate::egamma_HF;
03070 break;
03071 default:
03072 assert(0);
03073 }
03074
03075
03076 pfCandidates_->push_back( PFCandidate( charge,
03077 tmp,
03078 particleType ) );
03079
03080
03081
03082 pfCandidates_->back().
03083 setPositionAtECALEntrance(math::XYZPointF(cluster.position().X(),
03084 cluster.position().Y(),
03085 cluster.position().Z()));
03086
03087
03088 pfCandidates_->back().setVertex(vertexPos);
03089
03090 if(debug_)
03091 cout<<"** candidate: "<<pfCandidates_->back()<<endl;
03092
03093
03094 return pfCandidates_->size()-1;
03095
03096 }
03097
03098
03099
03100
03101 double
03102 PFAlgo::neutralHadronEnergyResolution(double clusterEnergyHCAL, double eta) const {
03103
03104
03105 if ( clusterEnergyHCAL < 1. ) clusterEnergyHCAL = 1.;
03106
03107 double resol = fabs(eta) < 1.48 ?
03108 sqrt (1.02*1.02/clusterEnergyHCAL + 0.065*0.065)
03109 :
03110 sqrt (1.20*1.20/clusterEnergyHCAL + 0.028*0.028);
03111
03112 return resol;
03113
03114 }
03115
03116 double
03117 PFAlgo::nSigmaHCAL(double clusterEnergyHCAL, double eta) const {
03118 double nS = fabs(eta) < 1.48 ?
03119 nSigmaHCAL_ * (1. + exp(-clusterEnergyHCAL/100.))
03120 :
03121 nSigmaHCAL_ * (1. + exp(-clusterEnergyHCAL/100.));
03122
03123 return nS;
03124 }
03125
03126
03127 ostream& operator<<(ostream& out, const PFAlgo& algo) {
03128 if(!out ) return out;
03129
03130 out<<"====== Particle Flow Algorithm ======= ";
03131 out<<endl;
03132 out<<"nSigmaECAL_ "<<algo.nSigmaECAL_<<endl;
03133 out<<"nSigmaHCAL_ "<<algo.nSigmaHCAL_<<endl;
03134 out<<endl;
03135 out<<(*(algo.calibration_))<<endl;
03136 out<<endl;
03137 out<<"reconstructed particles: "<<endl;
03138
03139 const std::auto_ptr< reco::PFCandidateCollection >&
03140 candidates = algo.pfCandidates();
03141
03142 if(!candidates.get() ) {
03143 out<<"candidates already transfered"<<endl;
03144 return out;
03145 }
03146 for(PFCandidateConstIterator ic=algo.pfCandidates_->begin();
03147 ic != algo.pfCandidates_->end(); ic++) {
03148 out<<(*ic)<<endl;
03149 }
03150
03151 return out;
03152 }
03153
03154 reco::PFBlockRef
03155 PFAlgo::createBlockRef( const reco::PFBlockCollection& blocks,
03156 unsigned bi ) {
03157
03158 if( blockHandle_.isValid() ) {
03159 return reco::PFBlockRef( blockHandle_, bi );
03160 }
03161 else {
03162 return reco::PFBlockRef( &blocks, bi );
03163 }
03164 }
03165
03166 void
03167 PFAlgo::associatePSClusters(unsigned iEcal,
03168 reco::PFBlockElement::Type psElementType,
03169 const reco::PFBlock& block,
03170 const edm::OwnVector< reco::PFBlockElement >& elements,
03171 const reco::PFBlock::LinkData& linkData,
03172 std::vector<bool>& active,
03173 std::vector<double>& psEne) {
03174
03175
03176
03177
03178
03179
03180
03181 std::multimap<double, unsigned> sortedPS;
03182 typedef std::multimap<double, unsigned>::iterator IE;
03183 block.associatedElements( iEcal, linkData,
03184 sortedPS, psElementType,
03185 reco::PFBlock::LINKTEST_ALL );
03186
03187
03188 double totalPS = 0.;
03189 for ( IE ips=sortedPS.begin(); ips!=sortedPS.end(); ++ips ) {
03190
03191
03192 unsigned iPS = ips->second;
03193
03194
03195
03196 if (!active[iPS]) continue;
03197
03198
03199 std::multimap<double, unsigned> sortedECAL;
03200 block.associatedElements( iPS, linkData,
03201 sortedECAL,
03202 reco::PFBlockElement::ECAL,
03203 reco::PFBlock::LINKTEST_ALL );
03204 unsigned jEcal = sortedECAL.begin()->second;
03205 if ( jEcal != iEcal ) continue;
03206
03207
03208 PFBlockElement::Type pstype = elements[ iPS ].type();
03209 assert( pstype == psElementType );
03210 PFClusterRef psclusterref = elements[iPS].clusterRef();
03211 assert(!psclusterref.isNull() );
03212 totalPS += psclusterref->energy();
03213 psEne[0] += psclusterref->energy();
03214 active[iPS] = false;
03215 }
03216
03217
03218 }
03219
03220
03221 bool
03222 PFAlgo::isFromSecInt(const reco::PFBlockElement& eTrack, string order) const {
03223
03224 reco::PFBlockElement::TrackType T_TO_DISP = reco::PFBlockElement::T_TO_DISP;
03225 reco::PFBlockElement::TrackType T_FROM_DISP = reco::PFBlockElement::T_FROM_DISP;
03226
03227 reco::PFBlockElement::TrackType T_FROM_V0 = reco::PFBlockElement::T_FROM_V0;
03228
03229 bool bPrimary = (order.find("primary") != string::npos);
03230 bool bSecondary = (order.find("secondary") != string::npos);
03231 bool bAll = (order.find("all") != string::npos);
03232
03233 bool isToDisp = usePFNuclearInteractions_ && eTrack.trackType(T_TO_DISP);
03234 bool isFromDisp = usePFNuclearInteractions_ && eTrack.trackType(T_FROM_DISP);
03235
03236 if (bPrimary && isToDisp) return true;
03237 if (bSecondary && isFromDisp ) return true;
03238 if (bAll && ( isToDisp || isFromDisp ) ) return true;
03239
03240
03241
03242
03243
03244 bool isFromDecay = (bAll || bSecondary) && usePFDecays_ && eTrack.trackType(T_FROM_V0);
03245
03246 return isFromDecay;
03247
03248
03249 }
03250
03251
03252 void
03253 PFAlgo::setEGElectronCollection(const reco::GsfElectronCollection & egelectrons) {
03254 if(useEGElectrons_ && pfele_) pfele_->setEGElectronCollection(egelectrons);
03255 }
03256
03257 void
03258 PFAlgo::postCleaning() {
03259
03260
03261 if ( !postHFCleaning_ ) return;
03262
03263
03264 double metX = 0.;
03265 double metY = 0.;
03266 double sumet = 0;
03267 std::vector<unsigned int> pfCandidatesToBeRemoved;
03268 for(unsigned i=0; i<pfCandidates_->size(); i++) {
03269 const PFCandidate& pfc = (*pfCandidates_)[i];
03270 metX += pfc.px();
03271 metY += pfc.py();
03272 sumet += pfc.pt();
03273 }
03274 double met2 = metX*metX+metY*metY;
03275
03276 double significance = std::sqrt(met2/sumet);
03277 double significanceCor = significance;
03278 if ( significance > minSignificance_ ) {
03279
03280 double metXCor = metX;
03281 double metYCor = metY;
03282 double sumetCor = sumet;
03283 double met2Cor = met2;
03284 double deltaPhi = 3.14159;
03285 double deltaPhiPt = 100.;
03286 bool next = true;
03287 unsigned iCor = 1E9;
03288
03289
03290 while ( next ) {
03291
03292 double metReduc = -1.;
03293
03294 for(unsigned i=0; i<pfCandidates_->size(); ++i) {
03295 const PFCandidate& pfc = (*pfCandidates_)[i];
03296
03297
03298 if ( pfc.particleId() != reco::PFCandidate::h_HF &&
03299 pfc.particleId() != reco::PFCandidate::egamma_HF ) continue;
03300
03301
03302 if ( pfc.pt() < minHFCleaningPt_ ) continue;
03303
03304
03305 bool skip = false;
03306 for(unsigned j=0; j<pfCandidatesToBeRemoved.size(); ++j) {
03307 if ( i == pfCandidatesToBeRemoved[j] ) skip = true;
03308 if ( skip ) break;
03309 }
03310 if ( skip ) continue;
03311
03312
03313 deltaPhi = std::acos((metX*pfc.px()+metY*pfc.py())/(pfc.pt()*std::sqrt(met2)));
03314 deltaPhiPt = deltaPhi*pfc.pt();
03315 if ( deltaPhiPt > maxDeltaPhiPt_ ) continue;
03316
03317
03318 double metXInt = metX - pfc.px();
03319 double metYInt = metY - pfc.py();
03320 double sumetInt = sumet - pfc.pt();
03321 double met2Int = metXInt*metXInt+metYInt*metYInt;
03322 if ( met2Int < met2Cor ) {
03323 metXCor = metXInt;
03324 metYCor = metYInt;
03325 metReduc = (met2-met2Int)/met2Int;
03326 met2Cor = met2Int;
03327 sumetCor = sumetInt;
03328 significanceCor = std::sqrt(met2Cor/sumetCor);
03329 iCor = i;
03330 }
03331 }
03332
03333
03334 if ( metReduc > minDeltaMet_ ) {
03335 pfCandidatesToBeRemoved.push_back(iCor);
03336 metX = metXCor;
03337 metY = metYCor;
03338 sumet = sumetCor;
03339 met2 = met2Cor;
03340 } else {
03341
03342 next = false;
03343 }
03344 }
03345
03346
03347 if ( significance - significanceCor > minSignificanceReduction_ &&
03348 significanceCor < maxSignificance_ ) {
03349 std::cout << "Significance reduction = " << significance << " -> "
03350 << significanceCor << " = " << significanceCor - significance
03351 << std::endl;
03352 for(unsigned j=0; j<pfCandidatesToBeRemoved.size(); ++j) {
03353 std::cout << "Removed : " << (*pfCandidates_)[pfCandidatesToBeRemoved[j]] << std::endl;
03354 pfCleanedCandidates_->push_back( (*pfCandidates_)[ pfCandidatesToBeRemoved[j] ] );
03355 (*pfCandidates_)[pfCandidatesToBeRemoved[j]].rescaleMomentum(1E-6);
03356
03357
03358 }
03359 }
03360 }
03361
03362 }
03363
03364 void
03365 PFAlgo::checkCleaning( const reco::PFRecHitCollection& cleanedHits ) {
03366
03367
03368
03369 if ( !cleanedHits.size() ) return;
03370
03371
03372 double metX = 0.;
03373 double metY = 0.;
03374 double sumet = 0;
03375 std::vector<unsigned int> hitsToBeAdded;
03376 for(unsigned i=0; i<pfCandidates_->size(); i++) {
03377 const PFCandidate& pfc = (*pfCandidates_)[i];
03378 metX += pfc.px();
03379 metY += pfc.py();
03380 sumet += pfc.pt();
03381 }
03382 double met2 = metX*metX+metY*metY;
03383 double met2_Original = met2;
03384
03385
03386
03387 double metXCor = metX;
03388 double metYCor = metY;
03389 double sumetCor = sumet;
03390 double met2Cor = met2;
03391 bool next = true;
03392 unsigned iCor = 1E9;
03393
03394
03395 while ( next ) {
03396
03397 double metReduc = -1.;
03398
03399 for(unsigned i=0; i<cleanedHits.size(); ++i) {
03400 const PFRecHit& hit = cleanedHits[i];
03401 double length = std::sqrt(hit.position().Mag2());
03402 double px = hit.energy() * hit.position().x()/length;
03403 double py = hit.energy() * hit.position().y()/length;
03404 double pt = std::sqrt(px*px + py*py);
03405
03406
03407 bool skip = false;
03408 for(unsigned j=0; j<hitsToBeAdded.size(); ++j) {
03409 if ( i == hitsToBeAdded[j] ) skip = true;
03410 if ( skip ) break;
03411 }
03412 if ( skip ) continue;
03413
03414
03415 double metXInt = metX + px;
03416 double metYInt = metY + py;
03417 double sumetInt = sumet + pt;
03418 double met2Int = metXInt*metXInt+metYInt*metYInt;
03419
03420
03421 if ( met2Int < met2Cor ) {
03422 metXCor = metXInt;
03423 metYCor = metYInt;
03424 metReduc = (met2-met2Int)/met2Int;
03425 met2Cor = met2Int;
03426 sumetCor = sumetInt;
03427
03428 iCor = i;
03429 }
03430 }
03431
03432
03433
03434 if ( metReduc > minDeltaMet_ ) {
03435 hitsToBeAdded.push_back(iCor);
03436 metX = metXCor;
03437 metY = metYCor;
03438 sumet = sumetCor;
03439 met2 = met2Cor;
03440 } else {
03441
03442 next = false;
03443 }
03444 }
03445
03446
03447 if ( std::sqrt(met2_Original) - std::sqrt(met2) > 5. ) {
03448 if ( debug_ ) {
03449 std::cout << hitsToBeAdded.size() << " hits were re-added " << std::endl;
03450 std::cout << "MET reduction = " << std::sqrt(met2_Original) << " -> "
03451 << std::sqrt(met2Cor) << " = " << std::sqrt(met2Cor) - std::sqrt(met2_Original)
03452 << std::endl;
03453 std::cout << "Added after cleaning check : " << std::endl;
03454 }
03455 for(unsigned j=0; j<hitsToBeAdded.size(); ++j) {
03456 const PFRecHit& hit = cleanedHits[hitsToBeAdded[j]];
03457 PFCluster cluster(hit.layer(), hit.energy(),
03458 hit.position().x(), hit.position().y(), hit.position().z() );
03459 reconstructCluster(cluster,hit.energy());
03460 if ( debug_ ) {
03461 std::cout << pfCandidates_->back() << ". time = " << hit.time() << std::endl;
03462 }
03463 }
03464 }
03465
03466 }
03467
03468
03469
03470 void PFAlgo::setElectronExtraRef(const edm::OrphanHandle<reco::PFCandidateElectronExtraCollection >& extrah) {
03471 if(!usePFElectrons_) return;
03472
03473 unsigned size=pfCandidates_->size();
03474
03475 for(unsigned ic=0;ic<size;++ic) {
03476
03477 if((*pfCandidates_)[ic].particleId()==PFCandidate::e) {
03478
03479 PFElectronExtraEqual myExtraEqual((*pfCandidates_)[ic].gsfTrackRef());
03480 std::vector<PFCandidateElectronExtra>::const_iterator it=find_if(pfElectronExtra_.begin(),pfElectronExtra_.end(),myExtraEqual);
03481 if(it!=pfElectronExtra_.end()) {
03482
03483 reco::PFCandidateElectronExtraRef theRef(extrah,it-pfElectronExtra_.begin());
03484 (*pfCandidates_)[ic].setPFElectronExtraRef(theRef);
03485 }
03486 else {
03487 (*pfCandidates_)[ic].setPFElectronExtraRef(PFCandidateElectronExtraRef());
03488 }
03489 }
03490 else
03491 {
03492 if((*pfCandidates_)[ic].trackRef().isNonnull()) {
03493 PFElectronExtraKfEqual myExtraEqual((*pfCandidates_)[ic].trackRef());
03494 std::vector<PFCandidateElectronExtra>::const_iterator it=find_if(pfElectronExtra_.begin(),pfElectronExtra_.end(),myExtraEqual);
03495 if(it!=pfElectronExtra_.end()) {
03496 (*pfCandidates_)[ic].set_mva_e_pi(it->mvaVariable(PFCandidateElectronExtra::MVA_MVA));
03497 reco::PFCandidateElectronExtraRef theRef(extrah,it-pfElectronExtra_.begin());
03498 (*pfCandidates_)[ic].setPFElectronExtraRef(theRef);
03499 (*pfCandidates_)[ic].setGsfTrackRef(it->gsfTrackRef());
03500 }
03501 }
03502 }
03503
03504 }
03505
03506 size=pfElectronCandidates_->size();
03507 for(unsigned ic=0;ic<size;++ic) {
03508
03509 if((*pfElectronCandidates_)[ic].particleId()==PFCandidate::e) {
03510
03511 PFElectronExtraEqual myExtraEqual((*pfElectronCandidates_)[ic].gsfTrackRef());
03512 std::vector<PFCandidateElectronExtra>::const_iterator it=find_if(pfElectronExtra_.begin(),pfElectronExtra_.end(),myExtraEqual);
03513 if(it!=pfElectronExtra_.end()) {
03514 reco::PFCandidateElectronExtraRef theRef(extrah,it-pfElectronExtra_.begin());
03515 (*pfElectronCandidates_)[ic].setPFElectronExtraRef(theRef);
03516
03517 }
03518 }
03519 }
03520
03521 }
03522 void PFAlgo::setPhotonExtraRef(const edm::OrphanHandle<reco::PFCandidatePhotonExtraCollection >& ph_extrah) {
03523 if(!usePFPhotons_) return;
03524 unsigned int size=pfCandidates_->size();
03525 unsigned int sizePhExtra = pfPhotonExtra_.size();
03526 for(unsigned ic=0;ic<size;++ic) {
03527
03528 if((*pfCandidates_)[ic].particleId()==PFCandidate::gamma && (*pfCandidates_)[ic].mva_nothing_gamma() > 0.99) {
03529 if((*pfCandidates_)[ic].superClusterRef().isNonnull()) {
03530 bool found = false;
03531 for(unsigned pcextra=0;pcextra<sizePhExtra;++pcextra) {
03532 if((*pfCandidates_)[ic].superClusterRef() == pfPhotonExtra_[pcextra].superClusterRef()) {
03533 reco::PFCandidatePhotonExtraRef theRef(ph_extrah,pcextra);
03534 (*pfCandidates_)[ic].setPFPhotonExtraRef(theRef);
03535 found = true;
03536 break;
03537 }
03538 }
03539 if(!found)
03540 (*pfCandidates_)[ic].setPFPhotonExtraRef((PFCandidatePhotonExtraRef()));
03541 }
03542 else {
03543 (*pfCandidates_)[ic].setPFPhotonExtraRef((PFCandidatePhotonExtraRef()));
03544 }
03545 }
03546 }
03547 }