File indexing completed on 2024-09-07 04:37:51
0001 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0002 #include "RecoParticleFlow/PFProducer/interface/PFAlgo.h"
0003 #include "RecoParticleFlow/PFProducer/interface/PFMuonAlgo.h"
0004 #include "RecoParticleFlow/PFProducer/interface/PFElectronExtraEqual.h"
0005 #include "RecoParticleFlow/PFTracking/interface/PFTrackAlgoTools.h"
0006
0007 #include "TDecompChol.h"
0008
0009 #include <numeric>
0010 #include <fstream>
0011
0012 using namespace std;
0013 using namespace reco;
0014
0015 PFAlgo::PFAlgo(double nSigmaECAL,
0016 double nSigmaHCAL,
0017 double nSigmaHFEM,
0018 double nSigmaHFHAD,
0019 std::vector<double> resolHF_square,
0020 PFEnergyCalibration& calibration,
0021 PFEnergyCalibrationHF& thepfEnergyCalibrationHF,
0022 const edm::ParameterSet& pset)
0023 : pfCandidates_(new PFCandidateCollection),
0024 nSigmaECAL_(nSigmaECAL),
0025 nSigmaHCAL_(nSigmaHCAL),
0026 nSigmaHFEM_(nSigmaHFEM),
0027 nSigmaHFHAD_(nSigmaHFHAD),
0028 resolHF_square_(resolHF_square),
0029 calibration_(calibration),
0030 thepfEnergyCalibrationHF_(thepfEnergyCalibrationHF),
0031 connector_() {
0032 const edm::ParameterSet pfMuonAlgoParams = pset.getParameter<edm::ParameterSet>("PFMuonAlgoParameters");
0033 bool postMuonCleaning = pset.getParameter<bool>("postMuonCleaning");
0034 pfmu_ = std::make_unique<PFMuonAlgo>(pfMuonAlgoParams, postMuonCleaning);
0035
0036
0037 assert(resolHF_square_.size() == 3);
0038
0039
0040 muonHCAL_ = pset.getParameter<std::vector<double>>("muon_HCAL");
0041 muonECAL_ = pset.getParameter<std::vector<double>>("muon_ECAL");
0042 muonHO_ = pset.getParameter<std::vector<double>>("muon_HO");
0043 assert(muonHCAL_.size() == 2 && muonECAL_.size() == 2 && muonHO_.size() == 2);
0044 nSigmaTRACK_ = pset.getParameter<double>("nsigma_TRACK");
0045 ptError_ = pset.getParameter<double>("pt_Error");
0046 factors45_ = pset.getParameter<std::vector<double>>("factors_45");
0047 assert(factors45_.size() == 2);
0048
0049
0050 goodTrackDeadHcal_ptErrRel_ = pset.getParameter<double>("goodTrackDeadHcal_ptErrRel");
0051 goodTrackDeadHcal_chi2n_ = pset.getParameter<double>("goodTrackDeadHcal_chi2n");
0052 goodTrackDeadHcal_layers_ = pset.getParameter<uint32_t>("goodTrackDeadHcal_layers");
0053 goodTrackDeadHcal_validFr_ = pset.getParameter<double>("goodTrackDeadHcal_validFr");
0054 goodTrackDeadHcal_dxy_ = pset.getParameter<double>("goodTrackDeadHcal_dxy");
0055
0056 goodPixelTrackDeadHcal_minEta_ = pset.getParameter<double>("goodPixelTrackDeadHcal_minEta");
0057 goodPixelTrackDeadHcal_maxPt_ = pset.getParameter<double>("goodPixelTrackDeadHcal_maxPt");
0058 goodPixelTrackDeadHcal_ptErrRel_ = pset.getParameter<double>("goodPixelTrackDeadHcal_ptErrRel");
0059 goodPixelTrackDeadHcal_chi2n_ = pset.getParameter<double>("goodPixelTrackDeadHcal_chi2n");
0060 goodPixelTrackDeadHcal_maxLost3Hit_ = pset.getParameter<int32_t>("goodPixelTrackDeadHcal_maxLost3Hit");
0061 goodPixelTrackDeadHcal_maxLost4Hit_ = pset.getParameter<int32_t>("goodPixelTrackDeadHcal_maxLost4Hit");
0062 goodPixelTrackDeadHcal_dxy_ = pset.getParameter<double>("goodPixelTrackDeadHcal_dxy");
0063 goodPixelTrackDeadHcal_dz_ = pset.getParameter<double>("goodPixelTrackDeadHcal_dz");
0064 }
0065
0066 PFMuonAlgo* PFAlgo::getPFMuonAlgo() { return pfmu_.get(); }
0067
0068 void PFAlgo::setEGammaParameters(bool use_EGammaFilters, bool useProtectionsForJetMET) {
0069 useEGammaFilters_ = use_EGammaFilters;
0070
0071 if (!useEGammaFilters_)
0072 return;
0073
0074 useProtectionsForJetMET_ = useProtectionsForJetMET;
0075 }
0076 void PFAlgo::setEGammaCollections(const edm::View<reco::PFCandidate>& pfEgammaCandidates,
0077 const edm::ValueMap<reco::GsfElectronRef>& valueMapGedElectrons,
0078 const edm::ValueMap<reco::PhotonRef>& valueMapGedPhotons) {
0079 if (useEGammaFilters_) {
0080 pfEgammaCandidates_ = &pfEgammaCandidates;
0081 valueMapGedElectrons_ = &valueMapGedElectrons;
0082 valueMapGedPhotons_ = &valueMapGedPhotons;
0083 }
0084 }
0085
0086 void PFAlgo::setPostHFCleaningParameters(bool postHFCleaning, const edm::ParameterSet& pfHFCleaningParams) {
0087 postHFCleaning_ = postHFCleaning;
0088 minHFCleaningPt_ = pfHFCleaningParams.getParameter<double>("minHFCleaningPt");
0089 minSignificance_ = pfHFCleaningParams.getParameter<double>("minSignificance");
0090 maxSignificance_ = pfHFCleaningParams.getParameter<double>("maxSignificance");
0091 minSignificanceReduction_ = pfHFCleaningParams.getParameter<double>("minSignificanceReduction");
0092 maxDeltaPhiPt_ = pfHFCleaningParams.getParameter<double>("maxDeltaPhiPt");
0093 minDeltaMet_ = pfHFCleaningParams.getParameter<double>("minDeltaMet");
0094 }
0095
0096 void PFAlgo::setDisplacedVerticesParameters(bool rejectTracks_Bad,
0097 bool rejectTracks_Step45,
0098 bool usePFNuclearInteractions,
0099 bool usePFConversions,
0100 bool usePFDecays,
0101 double dptRel_DispVtx) {
0102 rejectTracks_Bad_ = rejectTracks_Bad;
0103 rejectTracks_Step45_ = rejectTracks_Step45;
0104 usePFNuclearInteractions_ = usePFNuclearInteractions;
0105 usePFConversions_ = usePFConversions;
0106 usePFDecays_ = usePFDecays;
0107 dptRel_DispVtx_ = dptRel_DispVtx;
0108 }
0109
0110 void PFAlgo::setPFVertexParameters(bool useVertex, reco::VertexCollection const& primaryVertices) {
0111 useVertices_ = useVertex;
0112
0113
0114 pfmu_->setInputsForCleaning(primaryVertices);
0115
0116
0117 bool primaryVertexFound = false;
0118 nVtx_ = primaryVertices.size();
0119 for (auto const& vertex : primaryVertices) {
0120 if (vertex.isValid() && (!vertex.isFake())) {
0121 primaryVertex_ = vertex;
0122 primaryVertexFound = true;
0123 break;
0124 }
0125 }
0126
0127 useVertices_ = useVertex && primaryVertexFound;
0128 }
0129
0130 void PFAlgo::reconstructParticles(const reco::PFBlockHandle& blockHandle, PFEGammaFilters const* pfegamma) {
0131 auto const& blocks = *blockHandle;
0132
0133
0134 pfCandidates_->clear();
0135
0136 LogTrace("PFAlgo|reconstructParticles")
0137 << "start of function PFAlgo::reconstructParticles, blocks.size()=" << blocks.size();
0138
0139
0140 std::list<reco::PFBlockRef> hcalBlockRefs;
0141 std::list<reco::PFBlockRef> ecalBlockRefs;
0142 std::list<reco::PFBlockRef> hoBlockRefs;
0143 std::list<reco::PFBlockRef> otherBlockRefs;
0144
0145 for (unsigned i = 0; i < blocks.size(); ++i) {
0146 reco::PFBlockRef blockref = reco::PFBlockRef(blockHandle, i);
0147
0148 const reco::PFBlock& block = *blockref;
0149 const edm::OwnVector<reco::PFBlockElement>& elements = block.elements();
0150
0151 bool singleEcalOrHcal = false;
0152 if (elements.size() == 1) {
0153 if (elements[0].type() == reco::PFBlockElement::ECAL) {
0154 ecalBlockRefs.push_back(blockref);
0155 singleEcalOrHcal = true;
0156 }
0157 if (elements[0].type() == reco::PFBlockElement::HCAL) {
0158 hcalBlockRefs.push_back(blockref);
0159 singleEcalOrHcal = true;
0160 }
0161 if (elements[0].type() == reco::PFBlockElement::HO) {
0162
0163 hoBlockRefs.push_back(blockref);
0164 singleEcalOrHcal = true;
0165 }
0166 }
0167
0168 if (!singleEcalOrHcal) {
0169 otherBlockRefs.push_back(blockref);
0170 }
0171 }
0172
0173 LogTrace("PFAlgo|reconstructParticles")
0174 << "# Ecal blocks: " << ecalBlockRefs.size() << ", # Hcal blocks: " << hcalBlockRefs.size()
0175 << ", # HO blocks: " << hoBlockRefs.size() << ", # Other blocks: " << otherBlockRefs.size();
0176
0177
0178
0179
0180 unsigned nblcks = 0;
0181 for (auto const& other : otherBlockRefs) {
0182 LogTrace("PFAlgo|reconstructParticles") << "processBlock, Block number " << nblcks++;
0183 processBlock(other, hcalBlockRefs, ecalBlockRefs, pfegamma);
0184 }
0185
0186 std::list<reco::PFBlockRef> empty;
0187
0188 unsigned hblcks = 0;
0189
0190 for (auto const& hcal : hcalBlockRefs) {
0191 LogTrace("PFAlgo|reconstructParticles") << "processBlock, HCAL block number " << hblcks++;
0192 processBlock(hcal, empty, empty, pfegamma);
0193 }
0194
0195 unsigned eblcks = 0;
0196
0197 for (auto const& ecal : ecalBlockRefs) {
0198 LogTrace("PFAlgo|reconstructParticles") << "processBlock, ECAL block number " << eblcks++;
0199 processBlock(ecal, empty, empty, pfegamma);
0200 }
0201
0202
0203 pfCleanedCandidates_.clear();
0204
0205 if (postHFCleaning_) {
0206 postCleaning();
0207 }
0208
0209
0210 pfmu_->postClean(pfCandidates_.get());
0211
0212
0213 if (muonHandle_.isValid())
0214 pfmu_->addMissingMuons(muonHandle_, pfCandidates_.get());
0215
0216 LogTrace("PFAlgo|reconstructParticles")
0217 << "end of function PFAlgo::reconstructParticles, pfCandidates_->size()=" << pfCandidates_->size();
0218 }
0219
0220 void PFAlgo::egammaFilters(const reco::PFBlockRef& blockref,
0221 std::vector<bool>& active,
0222 PFEGammaFilters const* pfegamma) {
0223
0224
0225 unsigned int negmcandidates = pfEgammaCandidates_->size();
0226 LogTrace("PFAlgo|egammaFilters") << "start of function PFAlgo::egammaFilters(), negmcandidates=" << negmcandidates;
0227
0228 for (unsigned int ieg = 0; ieg < negmcandidates; ++ieg) {
0229
0230 reco::PFCandidateRef pfEgmRef = pfEgammaCandidates_->refAt(ieg).castTo<reco::PFCandidateRef>();
0231
0232 const PFCandidate::ElementsInBlocks& theElements = (*pfEgmRef).elementsInBlocks();
0233 PFCandidate::ElementsInBlocks::const_iterator iegfirst = theElements.begin();
0234 bool sameBlock = false;
0235 bool isGoodElectron = false;
0236 bool isGoodPhoton = false;
0237 bool isPrimaryElectron = false;
0238 if (iegfirst->first == blockref)
0239 sameBlock = true;
0240 if (sameBlock) {
0241 LogTrace("PFAlgo|egammaFilters") << " I am in looping on EGamma Candidates: pt " << (*pfEgmRef).pt()
0242 << " eta,phi " << (*pfEgmRef).eta() << ", " << (*pfEgmRef).phi() << " charge "
0243 << (*pfEgmRef).charge();
0244
0245 if ((*pfEgmRef).gsfTrackRef().isNonnull()) {
0246 reco::GsfElectronRef gedEleRef = (*valueMapGedElectrons_)[pfEgmRef];
0247 if (gedEleRef.isNonnull()) {
0248 isGoodElectron = pfegamma->passElectronSelection(*gedEleRef, *pfEgmRef, nVtx_);
0249 isPrimaryElectron = pfegamma->isElectron(*gedEleRef);
0250 if (isGoodElectron)
0251 LogTrace("PFAlgo|egammaFilters")
0252 << "** Good Electron, pt " << gedEleRef->pt() << " eta, phi " << gedEleRef->eta() << ", "
0253 << gedEleRef->phi() << " charge " << gedEleRef->charge() << " isPrimary " << isPrimaryElectron
0254 << " isoDr03 "
0255 << (gedEleRef->dr03TkSumPt() + gedEleRef->dr03EcalRecHitSumEt() + gedEleRef->dr03HcalTowerSumEt())
0256 << " mva_isolated " << gedEleRef->mva_Isolated() << " mva_e_pi " << gedEleRef->mva_e_pi();
0257 }
0258 }
0259 if ((*pfEgmRef).superClusterRef().isNonnull()) {
0260 reco::PhotonRef gedPhoRef = (*valueMapGedPhotons_)[pfEgmRef];
0261 if (gedPhoRef.isNonnull()) {
0262 isGoodPhoton = pfegamma->passPhotonSelection(*gedPhoRef);
0263 if (isGoodPhoton)
0264 LogTrace("PFAlgo|egammaFilters") << "** Good Photon, pt " << gedPhoRef->pt() << " eta, phi "
0265 << gedPhoRef->eta() << ", " << gedPhoRef->phi() << endl;
0266 }
0267 }
0268 }
0269
0270 if (isGoodElectron && isGoodPhoton) {
0271 if (isPrimaryElectron)
0272 isGoodPhoton = false;
0273 else
0274 isGoodElectron = false;
0275 }
0276
0277
0278 if (isGoodElectron) {
0279 reco::GsfElectronRef gedEleRef = (*valueMapGedElectrons_)[pfEgmRef];
0280 reco::PFCandidate myPFElectron = *pfEgmRef;
0281
0282 bool lockTracks = false;
0283 bool isSafe = true;
0284 if (useProtectionsForJetMET_) {
0285 lockTracks = true;
0286 isSafe = pfegamma->isElectronSafeForJetMET(*gedEleRef, myPFElectron, primaryVertex_, lockTracks);
0287 }
0288
0289 if (isSafe) {
0290 reco::PFCandidate::ParticleType particleType = reco::PFCandidate::e;
0291 myPFElectron.setParticleType(particleType);
0292 myPFElectron.setCharge(gedEleRef->charge());
0293 myPFElectron.setP4(gedEleRef->p4());
0294 myPFElectron.set_mva_e_pi(gedEleRef->mva_e_pi());
0295 myPFElectron.set_mva_Isolated(gedEleRef->mva_Isolated());
0296
0297 myPFElectron.set_dnn_e_sigIsolated(gedEleRef->dnn_signal_Isolated());
0298 myPFElectron.set_dnn_e_sigNonIsolated(gedEleRef->dnn_signal_nonIsolated());
0299 myPFElectron.set_dnn_e_bkgNonIsolated(gedEleRef->dnn_bkg_nonIsolated());
0300 myPFElectron.set_dnn_e_bkgTau(gedEleRef->dnn_bkg_Tau());
0301 myPFElectron.set_dnn_e_bkgPhoton(gedEleRef->dnn_bkg_Photon());
0302
0303 LogTrace("PFAlgo|egammaFilters") << " PFAlgo: found an electron with NEW EGamma code ";
0304 LogTrace("PFAlgo|egammaFilters") << " myPFElectron: pt " << myPFElectron.pt() << " eta,phi "
0305 << myPFElectron.eta() << ", " << myPFElectron.phi() << " mva "
0306 << myPFElectron.mva_e_pi() << " charge " << myPFElectron.charge();
0307
0308 LogTrace("PFAlgo|egammaFilters") << " THE BLOCK " << *blockref;
0309 for (auto const& eb : theElements) {
0310 active[eb.second] = false;
0311 LogTrace("PFAlgo|egammaFilters") << " Elements used " << eb.second;
0312 }
0313
0314
0315 if (lockTracks) {
0316 const PFCandidate::ElementsInBlocks& extraTracks = myPFElectron.egammaExtraRef()->extraNonConvTracks();
0317 for (auto const& trk : extraTracks) {
0318 LogTrace("PFAlgo|egammaFilters") << " Extra locked track " << trk.second;
0319 active[trk.second] = false;
0320 }
0321 }
0322
0323 LogTrace("PFAlgo|egammaFilters") << "Creating PF electron: pt=" << myPFElectron.pt()
0324 << " eta=" << myPFElectron.eta() << " phi=" << myPFElectron.phi();
0325 pfCandidates_->push_back(myPFElectron);
0326
0327 } else {
0328 LogTrace("PFAlgo|egammaFilters") << "PFAlgo: Electron DISCARDED, NOT SAFE FOR JETMET ";
0329 }
0330 }
0331
0332 if (isGoodPhoton) {
0333 reco::PhotonRef gedPhoRef = (*valueMapGedPhotons_)[pfEgmRef];
0334 reco::PFCandidate myPFPhoton = *pfEgmRef;
0335 bool isSafe = true;
0336 if (useProtectionsForJetMET_) {
0337 isSafe = pfegamma->isPhotonSafeForJetMET(*gedPhoRef, myPFPhoton);
0338 }
0339
0340 if (isSafe) {
0341 reco::PFCandidate::ParticleType particleType = reco::PFCandidate::gamma;
0342 myPFPhoton.setParticleType(particleType);
0343 myPFPhoton.setCharge(0);
0344 myPFPhoton.set_mva_nothing_gamma(1.);
0345 ::math::XYZPoint v(primaryVertex_.x(), primaryVertex_.y(), primaryVertex_.z());
0346 myPFPhoton.setVertex(v);
0347 myPFPhoton.setP4(gedPhoRef->p4());
0348
0349 myPFPhoton.set_dnn_gamma(gedPhoRef->pfDNN());
0350 LogTrace("PFAlgo|egammaFilters") << " PFAlgo: found a photon with NEW EGamma code ";
0351 LogTrace("PFAlgo|egammaFilters") << " myPFPhoton: pt " << myPFPhoton.pt() << " eta,phi " << myPFPhoton.eta()
0352 << ", " << myPFPhoton.phi() << " charge " << myPFPhoton.charge();
0353
0354
0355 LogTrace("PFAlgo|egammaFilters") << " THE BLOCK " << *blockref;
0356 for (auto const& eb : theElements) {
0357 active[eb.second] = false;
0358 LogTrace("PFAlgo|egammaFilters") << " Elements used " << eb.second;
0359 }
0360 LogTrace("PFAlgo|egammaFilters") << "Creating PF photon: pt=" << myPFPhoton.pt() << " eta=" << myPFPhoton.eta()
0361 << " phi=" << myPFPhoton.phi();
0362 pfCandidates_->push_back(myPFPhoton);
0363
0364 }
0365 }
0366 }
0367 LogTrace("PFAlgo|egammaFilters") << "end of function PFAlgo::egammaFilters";
0368 }
0369
0370 void PFAlgo::conversionAlgo(const edm::OwnVector<reco::PFBlockElement>& elements, std::vector<bool>& active) {
0371 LogTrace("PFAlgo|conversionAlgo") << "start of function PFAlgo::conversionAlgo(), elements.size()="
0372 << elements.size();
0373 for (unsigned iEle = 0; iEle < elements.size(); iEle++) {
0374 PFBlockElement::Type type = elements[iEle].type();
0375 if (type == PFBlockElement::TRACK) {
0376 LogTrace("PFAlgo|conversionAlgo") << "elements[" << iEle << "].type() == TRACK, active[" << iEle
0377 << "]=" << active[iEle];
0378 if (elements[iEle].trackRef()->algo() == reco::TrackBase::conversionStep) {
0379 active[iEle] = false;
0380 }
0381 if (elements[iEle].trackRef()->quality(reco::TrackBase::highPurity)) {
0382 LogTrace("PFAlgo|conversionAlgo") << "Track is high purity";
0383 continue;
0384 }
0385 const auto* trackRef = dynamic_cast<const reco::PFBlockElementTrack*>((&elements[iEle]));
0386 if (!(trackRef->trackType(reco::PFBlockElement::T_FROM_GAMMACONV))) {
0387 LogTrace("PFAlgo|conversionAlgo") << "!trackType(T_FROM_GAMMACONV)";
0388 continue;
0389 }
0390 if (!elements[iEle].convRefs().empty()) {
0391 active[iEle] = false;
0392 }
0393 LogTrace("PFAlgo|conversionAlgo") << "active[iEle=" << iEle << "]=" << active[iEle];
0394 }
0395 }
0396 LogTrace("PFAlgo|conversionAlgo") << "end of function PFAlgo::conversionAlgo";
0397 }
0398
0399 bool PFAlgo::recoTracksNotHCAL(const reco::PFBlock& block,
0400 reco::PFBlock::LinkData& linkData,
0401 const edm::OwnVector<reco::PFBlockElement>& elements,
0402 const reco::PFBlockRef& blockref,
0403 std::vector<bool>& active,
0404 bool goodTrackDeadHcal,
0405 bool hasDeadHcal,
0406 unsigned int iTrack,
0407 std::multimap<double, unsigned>& ecalElems,
0408 reco::TrackRef& trackRef) {
0409 LogTrace("PFAlgo|recoTracksNotHCAL") << "start of function PFAlgo::recoTracksNotHCAL(), now dealing with tracks "
0410 "linked to no HCAL clusters. Was HCal active? "
0411 << (!hasDeadHcal);
0412
0413
0414
0415
0416 double trackMomentum = elements[iTrack].trackRef()->p();
0417 LogTrace("PFAlgo|recoTracksNotHCAL") << elements[iTrack];
0418
0419
0420
0421
0422 bool thisIsAMuon = PFMuonAlgo::isMuon(elements[iTrack]);
0423 if (thisIsAMuon)
0424 trackMomentum = 0.;
0425
0426
0427
0428 if (!thisIsAMuon && elements[iTrack].trackRef()->ptError() > ptError_) {
0429 double deficit = trackMomentum;
0430 double resol = neutralHadronEnergyResolution(trackMomentum, elements[iTrack].trackRef()->eta());
0431 resol *= trackMomentum;
0432
0433 if (!ecalElems.empty()) {
0434 unsigned thisEcal = ecalElems.begin()->second;
0435 reco::PFClusterRef clusterRef = elements[thisEcal].clusterRef();
0436 bool useCluster = true;
0437 if (hasDeadHcal) {
0438 std::multimap<double, unsigned> sortedTracks;
0439 block.associatedElements(
0440 thisEcal, linkData, sortedTracks, reco::PFBlockElement::TRACK, reco::PFBlock::LINKTEST_ALL);
0441 useCluster = (sortedTracks.begin()->second == iTrack);
0442 }
0443 if (useCluster) {
0444 deficit -= clusterRef->energy();
0445 resol = neutralHadronEnergyResolution(trackMomentum, clusterRef->positionREP().Eta());
0446 resol *= trackMomentum;
0447 }
0448 }
0449
0450 bool isPrimary = isFromSecInt(elements[iTrack], "primary");
0451
0452 if (deficit > nSigmaTRACK_ * resol && !isPrimary && !goodTrackDeadHcal) {
0453 active[iTrack] = false;
0454 LogTrace("PFAlgo|recoTracksNotHCAL")
0455 << elements[iTrack] << std::endl
0456 << " deficit " << deficit << " > " << nSigmaTRACK_ << " x " << resol << " track pt " << trackRef->pt()
0457 << " +- " << trackRef->ptError() << " layers valid " << trackRef->hitPattern().trackerLayersWithMeasurement()
0458 << ", lost " << trackRef->hitPattern().trackerLayersWithoutMeasurement(HitPattern::TRACK_HITS)
0459 << ", lost outer " << trackRef->hitPattern().trackerLayersWithoutMeasurement(HitPattern::MISSING_OUTER_HITS)
0460 << ", lost inner " << trackRef->hitPattern().trackerLayersWithoutMeasurement(HitPattern::MISSING_INNER_HITS)
0461 << "(valid fraction " << trackRef->validFraction() << ")"
0462 << " chi2/ndf " << trackRef->normalizedChi2() << " |dxy| "
0463 << std::abs(trackRef->dxy(primaryVertex_.position())) << " +- " << trackRef->dxyError() << " |dz| "
0464 << std::abs(trackRef->dz(primaryVertex_.position())) << " +- " << trackRef->dzError()
0465 << "is probably a fake (1) --> lock the track";
0466 return true;
0467 }
0468 }
0469
0470
0471
0472
0473 std::vector<unsigned> tmpi;
0474 std::vector<unsigned> kTrack;
0475
0476
0477 double dpt = trackRef->ptError();
0478
0479 if (rejectTracks_Step45_ && ecalElems.empty() && trackMomentum > 30. && dpt > 0.5 &&
0480 (PFTrackAlgoTools::step45(trackRef->algo()) && !goodTrackDeadHcal)) {
0481 double dptRel = dpt / trackRef->pt() * 100;
0482 bool isPrimaryOrSecondary = isFromSecInt(elements[iTrack], "all");
0483
0484 if (isPrimaryOrSecondary && dptRel < dptRel_DispVtx_) {
0485 return true;
0486 };
0487 unsigned nHits = elements[iTrack].trackRef()->hitPattern().trackerLayersWithMeasurement();
0488 unsigned int NLostHit = trackRef->hitPattern().trackerLayersWithoutMeasurement(HitPattern::TRACK_HITS);
0489
0490 LogTrace("PFAlgo|recoTracksNotHCAL") << "A track (algo = " << trackRef->algo() << ") with momentum "
0491 << trackMomentum << " / " << elements[iTrack].trackRef()->pt() << " +/- "
0492 << dpt << " / " << elements[iTrack].trackRef()->eta()
0493 << " without any link to ECAL/HCAL and with " << nHits << " (" << NLostHit
0494 << ") hits (lost hits) has been cleaned";
0495
0496 active[iTrack] = false;
0497 return true;
0498 }
0499
0500 tmpi.push_back(reconstructTrack(elements[iTrack]));
0501
0502 kTrack.push_back(iTrack);
0503 active[iTrack] = false;
0504
0505
0506 if (ecalElems.empty()) {
0507 (*pfCandidates_)[tmpi[0]].setEcalEnergy(0., 0.);
0508 (*pfCandidates_)[tmpi[0]].setHcalEnergy(0., 0.);
0509 (*pfCandidates_)[tmpi[0]].setHoEnergy(0., 0.);
0510 (*pfCandidates_)[tmpi[0]].setPs1Energy(0);
0511 (*pfCandidates_)[tmpi[0]].setPs2Energy(0);
0512 (*pfCandidates_)[tmpi[0]].addElementInBlock(blockref, kTrack[0]);
0513 return true;
0514 }
0515
0516
0517 const unsigned int thisEcal = ecalElems.begin()->second;
0518 reco::PFClusterRef clusterRef = elements[thisEcal].clusterRef();
0519 LogTrace("PFAlgo|recoTracksNotHCAL") << " is associated to " << elements[thisEcal];
0520
0521
0522 if (thisIsAMuon) {
0523 (*pfCandidates_)[tmpi[0]].setEcalEnergy(clusterRef->energy(), std::min(clusterRef->energy(), muonECAL_[0]));
0524 (*pfCandidates_)[tmpi[0]].setHcalEnergy(0., 0.);
0525 (*pfCandidates_)[tmpi[0]].setHoEnergy(0., 0.);
0526 (*pfCandidates_)[tmpi[0]].setPs1Energy(0);
0527 (*pfCandidates_)[tmpi[0]].setPs2Energy(0);
0528 (*pfCandidates_)[tmpi[0]].addElementInBlock(blockref, kTrack[0]);
0529 }
0530
0531 double slopeEcal = 1.;
0532 bool connectedToEcal = false;
0533 unsigned iEcal = 99999;
0534 double calibEcal = 0.;
0535 double calibHcal = 0.;
0536 double totalEcal = thisIsAMuon ? -muonECAL_[0] : 0.;
0537
0538
0539 std::multimap<double, unsigned> sortedTracks;
0540 block.associatedElements(thisEcal, linkData, sortedTracks, reco::PFBlockElement::TRACK, reco::PFBlock::LINKTEST_ALL);
0541 LogTrace("PFAlgo|recoTracksNotHCAL") << "the closest ECAL cluster, id " << thisEcal << ", has " << sortedTracks.size()
0542 << " associated tracks, now processing them. ";
0543
0544 if (hasDeadHcal && !sortedTracks.empty()) {
0545
0546 LogTrace("PFAlgo|recoTracksNotHCAL") << " the closest track to ECAL " << thisEcal << " is "
0547 << sortedTracks.begin()->second << " (distance " << sortedTracks.begin()->first
0548 << ")";
0549 if (sortedTracks.begin()->second != iTrack) {
0550 LogTrace("PFAlgo|recoTracksNotHCAL")
0551 << " the closest track to ECAL " << thisEcal << " is " << sortedTracks.begin()->second
0552 << " which is not the one being processed. Will skip ECAL linking for this track";
0553 (*pfCandidates_)[tmpi[0]].setEcalEnergy(0., 0.);
0554 (*pfCandidates_)[tmpi[0]].setHcalEnergy(0., 0.);
0555 (*pfCandidates_)[tmpi[0]].setHoEnergy(0., 0.);
0556 (*pfCandidates_)[tmpi[0]].setPs1Energy(0);
0557 (*pfCandidates_)[tmpi[0]].setPs2Energy(0);
0558 (*pfCandidates_)[tmpi[0]].addElementInBlock(blockref, kTrack[0]);
0559 return true;
0560 } else {
0561 LogTrace("PFAlgo|recoTracksNotHCAL")
0562 << " the closest track to ECAL " << thisEcal << " is " << sortedTracks.begin()->second
0563 << " which is the one being processed. Will skip ECAL linking for all other track";
0564 sortedTracks.clear();
0565 }
0566 }
0567
0568 for (auto const& trk : sortedTracks) {
0569 unsigned jTrack = trk.second;
0570
0571
0572 if (!active[jTrack])
0573 continue;
0574
0575
0576 if (jTrack == iTrack)
0577 continue;
0578
0579
0580
0581
0582
0583
0584
0585 std::multimap<double, unsigned> sortedECAL;
0586 block.associatedElements(jTrack, linkData, sortedECAL, reco::PFBlockElement::ECAL, reco::PFBlock::LINKTEST_ALL);
0587 if (sortedECAL.begin()->second != thisEcal)
0588 continue;
0589
0590
0591 bool thatIsAMuon = PFMuonAlgo::isMuon(elements[jTrack]);
0592 LogTrace("PFAlgo|recoTracksNotHCAL") << " found track " << jTrack << (thatIsAMuon ? " (muon) " : " (non-muon)")
0593 << ", with distance = " << sortedECAL.begin()->first;
0594
0595
0596 bool rejectFake = false;
0597 reco::TrackRef trackRef = elements[jTrack].trackRef();
0598 if (!thatIsAMuon && trackRef->ptError() > ptError_) {
0599
0600
0601 double deficit = trackMomentum + trackRef->p() - clusterRef->energy();
0602 double resol =
0603 nSigmaTRACK_ * neutralHadronEnergyResolution(trackMomentum + trackRef->p(), clusterRef->positionREP().Eta());
0604 resol *= (trackMomentum + trackRef->p());
0605 if (deficit > nSigmaTRACK_ * resol && !goodTrackDeadHcal) {
0606 rejectFake = true;
0607 kTrack.push_back(jTrack);
0608 active[jTrack] = false;
0609 LogTrace("PFAlgo|recoTracksNotHCAL")
0610 << elements[jTrack] << std::endl
0611 << "is probably a fake (2) --> lock the track"
0612 << "(trackMomentum = " << trackMomentum << ", clusterEnergy = " << clusterRef->energy()
0613 << ", deficit = " << deficit << " > " << nSigmaTRACK_ << " x " << resol
0614 << " assuming neutralHadronEnergyResolution "
0615 << neutralHadronEnergyResolution(trackMomentum + trackRef->p(), clusterRef->positionREP().Eta()) << ") ";
0616 }
0617 }
0618 if (rejectFake)
0619 continue;
0620
0621
0622
0623
0624 if (!thatIsAMuon) {
0625 LogTrace("PFAlgo|recoTracksNotHCAL") << "Track momentum increased from " << trackMomentum << " GeV ";
0626 trackMomentum += trackRef->p();
0627 LogTrace("PFAlgo|recoTracksNotHCAL") << "to " << trackMomentum << " GeV.";
0628 LogTrace("PFAlgo|recoTracksNotHCAL") << "with " << elements[jTrack];
0629 } else {
0630 totalEcal -= muonECAL_[0];
0631 totalEcal = std::max(totalEcal, 0.);
0632 }
0633
0634
0635
0636 tmpi.push_back(reconstructTrack(elements[jTrack]));
0637
0638 kTrack.push_back(jTrack);
0639 active[jTrack] = false;
0640
0641 if (thatIsAMuon) {
0642 (*pfCandidates_)[tmpi.back()].setEcalEnergy(clusterRef->energy(), std::min(clusterRef->energy(), muonECAL_[0]));
0643 (*pfCandidates_)[tmpi.back()].setHcalEnergy(0., 0.);
0644 (*pfCandidates_)[tmpi.back()].setHoEnergy(0., 0.);
0645 (*pfCandidates_)[tmpi.back()].setPs1Energy(0);
0646 (*pfCandidates_)[tmpi.back()].setPs2Energy(0);
0647 (*pfCandidates_)[tmpi.back()].addElementInBlock(blockref, kTrack.back());
0648 }
0649 }
0650
0651 LogTrace("PFAlgo|recoTracksNotHCAL") << "Loop over all associated ECAL clusters";
0652
0653 for (auto const& ecal : ecalElems) {
0654 const unsigned index = ecal.second;
0655 const PFBlockElement::Type type = elements[index].type();
0656 assert(type == PFBlockElement::ECAL);
0657 LogTrace("PFAlgo|recoTracksNotHCAL") << elements[index];
0658
0659
0660 if (!active[index]) {
0661 LogTrace("PFAlgo|recoTracksNotHCAL") << "is not active - ignore ";
0662 continue;
0663 }
0664
0665
0666 block.associatedElements(index, linkData, sortedTracks, reco::PFBlockElement::TRACK, reco::PFBlock::LINKTEST_ALL);
0667
0668 const bool skip = std::none_of(
0669 kTrack.begin(), kTrack.end(), [&](unsigned iTrack) { return sortedTracks.begin()->second == iTrack; });
0670
0671 if (skip) {
0672 LogTrace("PFAlgo|recoTracksNotHCAL") << "is closer to another track - ignore ";
0673 continue;
0674 }
0675
0676
0677 const reco::PFClusterRef clusterRef = elements[index].clusterRef();
0678 assert(!clusterRef.isNull());
0679
0680 LogTrace("PFAlgo|recoTracksNotHCAL") << "Ecal cluster with raw energy = " << clusterRef->energy()
0681 << " linked with distance = " << ecal.first;
0682
0683
0684
0685
0686
0687 vector<double> ps1Ene{0.};
0688 associatePSClusters(index, reco::PFBlockElement::PS1, block, elements, linkData, active, ps1Ene);
0689 vector<double> ps2Ene{0.};
0690 associatePSClusters(index, reco::PFBlockElement::PS2, block, elements, linkData, active, ps2Ene);
0691
0692
0693 const double ecalEnergy = clusterRef->energy();
0694 const double ecalEnergyCalibrated = clusterRef->correctedEnergy();
0695 LogTrace("PFAlgo|recoTracksNotHCAL") << "Corrected ECAL(+PS) energy = " << ecalEnergy;
0696
0697
0698
0699 totalEcal += ecalEnergy;
0700 const double previousCalibEcal = calibEcal;
0701 const double previousSlopeEcal = slopeEcal;
0702 calibEcal = std::max(totalEcal, 0.);
0703 calibHcal = 0.;
0704 calibration_.energyEmHad(
0705 trackMomentum, calibEcal, calibHcal, clusterRef->positionREP().Eta(), clusterRef->positionREP().Phi());
0706 if (totalEcal > 0.)
0707 slopeEcal = calibEcal / totalEcal;
0708
0709 LogTrace("PFAlgo|recoTracksNotHCAL") << "The total calibrated energy so far amounts to = " << calibEcal
0710 << " (slope = " << slopeEcal << ")";
0711
0712
0713 if (connectedToEcal && calibEcal - trackMomentum >= 0.) {
0714
0715
0716 calibEcal = previousCalibEcal;
0717 slopeEcal = previousSlopeEcal;
0718 totalEcal = calibEcal / slopeEcal;
0719
0720
0721
0722 active[index] = false;
0723
0724
0725 std::multimap<double, unsigned> assTracks;
0726 block.associatedElements(index, linkData, assTracks, reco::PFBlockElement::TRACK, reco::PFBlock::LINKTEST_ALL);
0727
0728 auto& ecalCand = (*pfCandidates_)[reconstructCluster(
0729 *clusterRef, ecalEnergyCalibrated)];
0730 ecalCand.setEcalEnergy(clusterRef->energy(), ecalEnergyCalibrated);
0731 ecalCand.setHcalEnergy(0., 0.);
0732 ecalCand.setHoEnergy(0., 0.);
0733 ecalCand.setPs1Energy(ps1Ene[0]);
0734 ecalCand.setPs2Energy(ps2Ene[0]);
0735 ecalCand.addElementInBlock(blockref, index);
0736
0737 if (!assTracks.empty()) {
0738 ecalCand.addElementInBlock(blockref, assTracks.begin()->second);
0739
0740
0741 const ::math::XYZPointF& chargedPosition =
0742 dynamic_cast<const reco::PFBlockElementTrack*>(&elements[assTracks.begin()->second])
0743 ->positionAtECALEntrance();
0744 ecalCand.setPositionAtECALEntrance(chargedPosition);
0745 }
0746 break;
0747 }
0748
0749
0750 connectedToEcal = true;
0751 iEcal = index;
0752 active[index] = false;
0753 for (unsigned ic : tmpi)
0754 (*pfCandidates_)[ic].addElementInBlock(blockref, iEcal);
0755
0756 }
0757
0758 bool bNeutralProduced = false;
0759
0760
0761 if (connectedToEcal) {
0762 reco::PFClusterRef pivotalRef = elements[iEcal].clusterRef();
0763
0764 double neutralEnergy = calibEcal - trackMomentum;
0765
0766
0767
0768
0769
0770
0771
0772
0773
0774
0775
0776
0777
0778
0779
0780
0781
0782
0783
0784
0785
0786
0787
0788
0789
0790
0791
0792
0793
0794
0795
0796
0797
0798
0799
0800 double resol = neutralHadronEnergyResolution(trackMomentum, pivotalRef->positionREP().Eta());
0801 resol *= trackMomentum;
0802 if (neutralEnergy > std::max(0.5, nSigmaECAL_ * resol)) {
0803 neutralEnergy /= slopeEcal;
0804 unsigned tmpj = reconstructCluster(*pivotalRef, neutralEnergy);
0805 (*pfCandidates_)[tmpj].setEcalEnergy(pivotalRef->energy(), neutralEnergy);
0806 (*pfCandidates_)[tmpj].setHcalEnergy(0., 0.);
0807 (*pfCandidates_)[tmpj].setHoEnergy(0., 0.);
0808 (*pfCandidates_)[tmpj].setPs1Energy(0.);
0809 (*pfCandidates_)[tmpj].setPs2Energy(0.);
0810 (*pfCandidates_)[tmpj].addElementInBlock(blockref, iEcal);
0811 bNeutralProduced = true;
0812 for (unsigned ic = 0; ic < kTrack.size(); ++ic)
0813 (*pfCandidates_)[tmpj].addElementInBlock(blockref, kTrack[ic]);
0814 }
0815
0816
0817 for (unsigned ic = 0; ic < tmpi.size(); ++ic) {
0818
0819 if ((*pfCandidates_)[tmpi[ic]].particleId() == reco::PFCandidate::mu)
0820 continue;
0821
0822 double fraction = trackMomentum > 0 ? (*pfCandidates_)[tmpi[ic]].trackRef()->p() / trackMomentum : 0;
0823 double ecalCal = bNeutralProduced ? (calibEcal - neutralEnergy * slopeEcal) * fraction : calibEcal * fraction;
0824 double ecalRaw = totalEcal * fraction;
0825
0826 LogTrace("PFAlgo|recoTracksNotHCAL")
0827 << "The fraction after photon supression is " << fraction << " calibrated ecal = " << ecalCal;
0828
0829 (*pfCandidates_)[tmpi[ic]].setEcalEnergy(ecalRaw, ecalCal);
0830 (*pfCandidates_)[tmpi[ic]].setHcalEnergy(0., 0.);
0831 (*pfCandidates_)[tmpi[ic]].setHoEnergy(0., 0.);
0832 (*pfCandidates_)[tmpi[ic]].setPs1Energy(0);
0833 (*pfCandidates_)[tmpi[ic]].setPs2Energy(0);
0834 (*pfCandidates_)[tmpi[ic]].addElementInBlock(blockref, kTrack[ic]);
0835 }
0836
0837 }
0838
0839
0840 for (unsigned ic = 0; ic < tmpi.size(); ++ic) {
0841 const PFCandidate& pfc = (*pfCandidates_)[tmpi[ic]];
0842 const PFCandidate::ElementsInBlocks& eleInBlocks = pfc.elementsInBlocks();
0843 if (eleInBlocks.empty()) {
0844 LogTrace("PFAlgo|recoTracksNotHCAL") << "Single track / Fill element in block! ";
0845 (*pfCandidates_)[tmpi[ic]].addElementInBlock(blockref, kTrack[ic]);
0846 }
0847 }
0848 LogTrace("PFAlgo|recoTracksNotHCAL") << "end of function PFAlgo::recoTracksNotHCAL";
0849 return false;
0850 }
0851
0852
0853
0854
0855 bool PFAlgo::checkAndReconstructSecondaryInteraction(const reco::PFBlockRef& blockref,
0856 const edm::OwnVector<reco::PFBlockElement>& elements,
0857 bool isActive,
0858 int iElement) {
0859 bool ret = isActive;
0860 if (isActive && isFromSecInt(elements[iElement], "primary")) {
0861 bool isPrimaryTrack =
0862 elements[iElement].displacedVertexRef(PFBlockElement::T_TO_DISP)->displacedVertexRef()->isTherePrimaryTracks();
0863 if (isPrimaryTrack) {
0864 LogTrace("PFAlgo|elementLoop") << "Primary Track reconstructed alone";
0865
0866 unsigned tmpi = reconstructTrack(elements[iElement]);
0867 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iElement);
0868 ret = false;
0869 }
0870 }
0871
0872 return ret;
0873 }
0874
0875 bool PFAlgo::checkHasDeadHcal(const std::multimap<double, unsigned>& hcalElems, const std::vector<bool>& deadArea) {
0876
0877
0878
0879
0880
0881
0882
0883
0884
0885
0886
0887
0888 bool hasDeadHcal = false;
0889 if (!hcalElems.empty() && deadArea[hcalElems.begin()->second]) {
0890 hasDeadHcal = true;
0891 }
0892
0893
0894
0895
0896
0897
0898 return hasDeadHcal;
0899 }
0900
0901
0902 bool PFAlgo::checkGoodTrackDeadHcal(const reco::TrackRef& trackRef, bool hasDeadHcal) {
0903 bool goodTrackDeadHcal = false;
0904 if (hasDeadHcal) {
0905 goodTrackDeadHcal = (trackRef->ptError() < goodTrackDeadHcal_ptErrRel_ * trackRef->pt() &&
0906 trackRef->normalizedChi2() < goodTrackDeadHcal_chi2n_ &&
0907 trackRef->hitPattern().trackerLayersWithMeasurement() >= goodTrackDeadHcal_layers_ &&
0908 trackRef->validFraction() > goodTrackDeadHcal_validFr_ &&
0909 std::abs(trackRef->dxy(primaryVertex_.position())) < goodTrackDeadHcal_dxy_);
0910
0911 if (!goodTrackDeadHcal && std::abs(trackRef->eta()) > goodPixelTrackDeadHcal_minEta_ &&
0912 trackRef->hitPattern().pixelLayersWithMeasurement() >= 3 &&
0913 trackRef->hitPattern().trackerLayersWithoutMeasurement(HitPattern::TRACK_HITS) == 0 &&
0914 trackRef->hitPattern().trackerLayersWithoutMeasurement(HitPattern::MISSING_INNER_HITS) == 0 &&
0915 trackRef->hitPattern().trackerLayersWithoutMeasurement(HitPattern::MISSING_OUTER_HITS) <=
0916 (trackRef->hitPattern().pixelLayersWithMeasurement() > 3 ? goodPixelTrackDeadHcal_maxLost4Hit_
0917 : goodPixelTrackDeadHcal_maxLost3Hit_) &&
0918 trackRef->normalizedChi2() < goodPixelTrackDeadHcal_chi2n_ &&
0919 std::abs(trackRef->dxy(primaryVertex_.position())) < goodPixelTrackDeadHcal_dxy_ &&
0920 std::abs(trackRef->dz(primaryVertex_.position())) < goodPixelTrackDeadHcal_dz_ &&
0921 trackRef->ptError() < goodPixelTrackDeadHcal_ptErrRel_ * trackRef->pt() &&
0922 trackRef->pt() < goodPixelTrackDeadHcal_maxPt_) {
0923 goodTrackDeadHcal = true;
0924
0925 }
0926
0927 LogTrace("PFAlgo|elementLoop")
0928 << " track pt " << trackRef->pt() << " +- " << trackRef->ptError() << " layers valid "
0929 << trackRef->hitPattern().trackerLayersWithMeasurement() << ", lost "
0930 << trackRef->hitPattern().trackerLayersWithoutMeasurement(HitPattern::TRACK_HITS) << ", lost outer "
0931 << trackRef->hitPattern().trackerLayersWithoutMeasurement(HitPattern::MISSING_OUTER_HITS) << ", lost inner "
0932 << trackRef->hitPattern().trackerLayersWithoutMeasurement(HitPattern::MISSING_INNER_HITS) << "(valid fraction "
0933 << trackRef->validFraction() << ")"
0934 << " chi2/ndf " << trackRef->normalizedChi2() << " |dxy| " << std::abs(trackRef->dxy(primaryVertex_.position()))
0935 << " +- " << trackRef->dxyError() << " |dz| " << std::abs(trackRef->dz(primaryVertex_.position())) << " +- "
0936 << trackRef->dzError() << (goodTrackDeadHcal ? " passes " : " fails ") << "quality cuts";
0937 }
0938 return goodTrackDeadHcal;
0939 }
0940
0941 void PFAlgo::relinkTrackToHcal(const reco::PFBlock& block,
0942 std::multimap<double, unsigned>& ecalElems,
0943 std::multimap<double, unsigned>& hcalElems,
0944 const std::vector<bool>& active,
0945 reco::PFBlock::LinkData& linkData,
0946 unsigned int iTrack) {
0947 unsigned index = ecalElems.begin()->second;
0948 std::multimap<double, unsigned> sortedTracks;
0949 block.associatedElements(index, linkData, sortedTracks, reco::PFBlockElement::TRACK, reco::PFBlock::LINKTEST_ALL);
0950 LogTrace("PFAlgo|elementLoop") << "The closest ECAL cluster is linked to " << sortedTracks.size()
0951 << " tracks, with distance = " << ecalElems.begin()->first;
0952
0953 LogTrace("PFAlgo|elementLoop") << "Looping over sortedTracks";
0954
0955 for (auto const& trk : sortedTracks) {
0956 unsigned jTrack = trk.second;
0957 LogTrace("PFAlgo|elementLoop") << "jTrack=" << jTrack;
0958
0959 if (!active[jTrack])
0960 continue;
0961 LogTrace("PFAlgo|elementLoop") << "active[jTrack]=" << active[jTrack];
0962
0963
0964 if (jTrack == iTrack)
0965 continue;
0966 LogTrace("PFAlgo|elementLoop") << "skipping jTrack=" << jTrack << " for same iTrack";
0967
0968
0969
0970 std::multimap<double, unsigned> sortedECAL;
0971 block.associatedElements(jTrack, linkData, sortedECAL, reco::PFBlockElement::ECAL, reco::PFBlock::LINKTEST_ALL);
0972 if (sortedECAL.begin()->second != index)
0973 continue;
0974 LogTrace("PFAlgo|elementLoop") << " track " << jTrack << " with closest ECAL identical ";
0975
0976
0977 std::multimap<double, unsigned> sortedHCAL;
0978 block.associatedElements(jTrack, linkData, sortedHCAL, reco::PFBlockElement::HCAL, reco::PFBlock::LINKTEST_ALL);
0979 if (sortedHCAL.empty())
0980 continue;
0981 LogTrace("PFAlgo|elementLoop") << " and with an HCAL cluster " << sortedHCAL.begin()->second;
0982
0983
0984 block.setLink(iTrack, sortedHCAL.begin()->second, sortedECAL.begin()->first, linkData, PFBlock::LINKTEST_RECHIT);
0985
0986 }
0987
0988
0989 block.associatedElements(iTrack, linkData, hcalElems, reco::PFBlockElement::HCAL, reco::PFBlock::LINKTEST_ALL);
0990
0991 if (!hcalElems.empty())
0992 LogTrace("PFAlgo|elementLoop") << "Track linked back to HCAL due to ECAL sharing with other tracks";
0993 }
0994
0995 void PFAlgo::elementLoop(const reco::PFBlock& block,
0996 reco::PFBlock::LinkData& linkData,
0997 const edm::OwnVector<reco::PFBlockElement>& elements,
0998 std::vector<bool>& active,
0999 const reco::PFBlockRef& blockref,
1000 ElementIndices& inds,
1001 std::vector<bool>& deadArea) {
1002 LogTrace("PFAlgo|elementLoop") << "start of function PFAlgo::elementLoop, elements.size()" << elements.size();
1003
1004 for (unsigned iEle = 0; iEle < elements.size(); iEle++) {
1005 PFBlockElement::Type type = elements[iEle].type();
1006
1007 LogTrace("PFAlgo|elementLoop") << "elements[iEle=" << iEle << "]=" << elements[iEle];
1008
1009
1010 auto ret_decideType = decideType(elements, type, active, inds, deadArea, iEle);
1011 if (ret_decideType == 1) {
1012 LogTrace("PFAlgo|elementLoop") << "ret_decideType==1, continuing";
1013 continue;
1014 }
1015 LogTrace("PFAlgo|elementLoop") << "ret_decideType=" << ret_decideType << " type=" << type;
1016
1017 active[iEle] = checkAndReconstructSecondaryInteraction(blockref, elements, active[iEle], iEle);
1018
1019 if (!active[iEle]) {
1020 LogTrace("PFAlgo|elementLoop") << "Already used by electrons, muons, conversions";
1021 continue;
1022 }
1023
1024 reco::TrackRef trackRef = elements[iEle].trackRef();
1025 assert(!trackRef.isNull());
1026
1027 LogTrace("PFAlgo|elementLoop") << "PFAlgo:processBlock"
1028 << " trackIs.size()=" << inds.trackIs.size()
1029 << " ecalIs.size()=" << inds.ecalIs.size() << " hcalIs.size()=" << inds.hcalIs.size()
1030 << " hoIs.size()=" << inds.hoIs.size() << " hfEmIs.size()=" << inds.hfEmIs.size()
1031 << " hfHadIs.size()=" << inds.hfHadIs.size();
1032
1033
1034
1035
1036
1037 std::multimap<double, unsigned> ecalElems;
1038 block.associatedElements(iEle, linkData, ecalElems, reco::PFBlockElement::ECAL, reco::PFBlock::LINKTEST_ALL);
1039
1040 std::multimap<double, unsigned> hcalElems;
1041 block.associatedElements(iEle, linkData, hcalElems, reco::PFBlockElement::HCAL, reco::PFBlock::LINKTEST_ALL);
1042
1043 std::multimap<double, unsigned> hfEmElems;
1044 std::multimap<double, unsigned> hfHadElems;
1045 block.associatedElements(iEle, linkData, hfEmElems, reco::PFBlockElement::HFEM, reco::PFBlock::LINKTEST_ALL);
1046 block.associatedElements(iEle, linkData, hfHadElems, reco::PFBlockElement::HFHAD, reco::PFBlock::LINKTEST_ALL);
1047
1048 LogTrace("PFAlgo|elementLoop") << "\tTrack " << iEle << " is linked to " << ecalElems.size() << " ecal and "
1049 << hcalElems.size() << " hcal and " << hfEmElems.size() << " hfEm and "
1050 << hfHadElems.size() << " hfHad elements";
1051
1052 #ifdef EDM_ML_DEBUG
1053 for (const auto& pair : ecalElems) {
1054 LogTrace("PFAlgo|elementLoop") << "ecal: dist " << pair.first << "\t elem " << pair.second;
1055 }
1056 for (const auto& pair : hcalElems) {
1057 LogTrace("PFAlgo|elementLoop") << "hcal: dist " << pair.first << "\t elem " << pair.second
1058 << (deadArea[pair.second] ? " DEAD AREA MARKER" : "");
1059 }
1060 #endif
1061
1062 const bool hasDeadHcal = checkHasDeadHcal(hcalElems, deadArea);
1063 if (hasDeadHcal) {
1064 hcalElems.clear();
1065 }
1066 const bool goodTrackDeadHcal = checkGoodTrackDeadHcal(trackRef, hasDeadHcal);
1067
1068
1069
1070
1071 if (hcalElems.empty() && !ecalElems.empty() && !hasDeadHcal) {
1072 relinkTrackToHcal(block, ecalElems, hcalElems, active, linkData, iEle);
1073 }
1074
1075
1076
1077
1078
1079
1080
1081
1082 std::multimap<double, unsigned> gsfElems;
1083 block.associatedElements(iEle, linkData, gsfElems, reco::PFBlockElement::GSF);
1084
1085 if (hcalElems.empty()) {
1086 LogTrace("PFAlgo|elementLoop") << "no hcal element connected to track " << iEle;
1087 }
1088
1089
1090 bool hcalFound = false;
1091 bool hfhadFound = false;
1092
1093 LogTrace("PFAlgo|elementLoop") << "now looping on elements associated to the track: ecalElems";
1094
1095
1096
1097 for (auto const& ecal : ecalElems) {
1098 unsigned index = ecal.second;
1099
1100 PFBlockElement::Type type = elements[index].type();
1101 #ifdef EDM_ML_DEBUG
1102 double dist = ecal.first;
1103 LogTrace("PFAlgo|elementLoop") << "\telement " << elements[index] << " linked with distance = " << dist;
1104 if (!active[index])
1105 LogTrace("PFAlgo|elementLoop") << "This ECAL is already used - skip it";
1106 #endif
1107 assert(type == PFBlockElement::ECAL);
1108
1109
1110 if (!active[index])
1111 continue;
1112
1113
1114
1115
1116
1117 if (!hcalElems.empty())
1118 LogTrace("PFAlgo|elementLoop") << "\t\tat least one hcal element connected to the track."
1119 << " Sparing Ecal cluster for the hcal loop";
1120
1121 }
1122
1123
1124
1125
1126 if (hcalElems.empty() && hfHadElems.empty()) {
1127 auto ret_continue = recoTracksNotHCAL(
1128 block, linkData, elements, blockref, active, goodTrackDeadHcal, hasDeadHcal, iEle, ecalElems, trackRef);
1129 if (ret_continue) {
1130 continue;
1131 }
1132 }
1133
1134
1135
1136 for (auto const& hcal : hcalElems) {
1137 unsigned index = hcal.second;
1138
1139 PFBlockElement::Type type = elements[index].type();
1140
1141 #ifdef EDM_ML_DEBUG
1142 double dist = block.dist(iEle, index, linkData, reco::PFBlock::LINKTEST_ALL);
1143 LogTrace("PFAlgo|elementLoop") << "\telement " << elements[index] << " linked with distance " << dist;
1144 #endif
1145 assert(type == PFBlockElement::HCAL);
1146
1147
1148
1149 if (!hcalFound) {
1150 LogTrace("PFAlgo|elementLoop") << "\t\tclosest hcal cluster, doing nothing";
1151
1152 hcalFound = true;
1153
1154
1155
1156 } else {
1157
1158 LogTrace("PFAlgo|elementLoop") << "\t\tsecondary hcal cluster. unlinking";
1159 block.setLink(iEle, index, -1., linkData, PFBlock::LINKTEST_RECHIT);
1160 }
1161 }
1162
1163
1164
1165
1166 for (auto const& hfhad : hfHadElems) {
1167 unsigned index = hfhad.second;
1168
1169 PFBlockElement::Type type = elements[index].type();
1170
1171 #ifdef EDM_ML_DEBUG
1172 double dist = block.dist(iEle, index, linkData, reco::PFBlock::LINKTEST_ALL);
1173 LogTrace("PFAlgo|elementLoop") << "\telement " << elements[index] << " linked with distance " << dist;
1174 #endif
1175 assert(type == PFBlockElement::HFHAD);
1176
1177
1178
1179 if (!hfhadFound) {
1180 LogTrace("PFAlgo|elementLoop") << "\t\tclosest hfhad cluster, doing nothing";
1181
1182 hfhadFound = true;
1183
1184 } else {
1185
1186 LogTrace("PFAlgo|elementLoop") << "\t\tsecondary hfhad cluster. unlinking";
1187 block.setLink(iEle, index, -1., linkData, PFBlock::LINKTEST_RECHIT);
1188 }
1189 }
1190
1191 LogTrace("PFAlgo|elementLoop") << "end of loop over iEle";
1192 }
1193 LogTrace("PFAlgo|elementLoop") << "end of function PFAlgo::elementLoop";
1194 }
1195
1196
1197
1198
1199 int PFAlgo::decideType(const edm::OwnVector<reco::PFBlockElement>& elements,
1200 const reco::PFBlockElement::Type type,
1201 std::vector<bool>& active,
1202 ElementIndices& inds,
1203 std::vector<bool>& deadArea,
1204 unsigned int iEle) {
1205 switch (type) {
1206 case PFBlockElement::TRACK:
1207 if (active[iEle]) {
1208 inds.trackIs.push_back(iEle);
1209 LogTrace("PFAlgo|decideType") << "TRACK, stored index, continue";
1210 }
1211 break;
1212 case PFBlockElement::ECAL:
1213 if (active[iEle]) {
1214 inds.ecalIs.push_back(iEle);
1215 LogTrace("PFAlgo|decideType") << "ECAL, stored index, continue";
1216 }
1217 return 1;
1218 case PFBlockElement::HCAL:
1219 if (active[iEle]) {
1220 if (elements[iEle].clusterRef()->flags() & reco::CaloCluster::badHcalMarker) {
1221 LogTrace("PFAlgo|decideType") << "HCAL DEAD AREA: remember and skip.";
1222 active[iEle] = false;
1223 deadArea[iEle] = true;
1224 return 1;
1225 }
1226 inds.hcalIs.push_back(iEle);
1227 LogTrace("PFAlgo|decideType") << "HCAL, stored index, continue";
1228 }
1229 return 1;
1230 case PFBlockElement::HO:
1231 if (useHO_) {
1232 if (active[iEle]) {
1233 inds.hoIs.push_back(iEle);
1234 LogTrace("PFAlgo|decideType") << "HO, stored index, continue";
1235 }
1236 }
1237 return 1;
1238 case PFBlockElement::HFEM:
1239 if (active[iEle]) {
1240 inds.hfEmIs.push_back(iEle);
1241 LogTrace("PFAlgo|decideType") << "HFEM, stored index, continue";
1242 }
1243 return 1;
1244 case PFBlockElement::HFHAD:
1245 if (active[iEle]) {
1246 inds.hfHadIs.push_back(iEle);
1247 LogTrace("PFAlgo|decideType") << "HFHAD, stored index, continue";
1248 }
1249 return 1;
1250 default:
1251 return 1;
1252 }
1253 LogTrace("PFAlgo|decideType") << "Did not match type to anything, return 0";
1254 return 0;
1255 }
1256
1257 void PFAlgo::createCandidatesHF(const reco::PFBlock& block,
1258 reco::PFBlock::LinkData& linkData,
1259 const edm::OwnVector<reco::PFBlockElement>& elements,
1260 std::vector<bool>& active,
1261 const reco::PFBlockRef& blockref,
1262 ElementIndices& inds) {
1263 LogTrace("PFAlgo|createCandidatesHF") << "starting function PFAlgo::createCandidatesHF";
1264
1265 bool trackInBlock = !inds.trackIs.empty();
1266
1267
1268 if (!trackInBlock)
1269 for (unsigned iEle = 0; iEle < elements.size(); iEle++) {
1270 PFBlockElement::Type type = elements[iEle].type();
1271 if (type == PFBlockElement::TRACK) {
1272 trackInBlock = true;
1273 break;
1274 }
1275 }
1276
1277
1278 if (!trackInBlock)
1279 assert(inds.hfEmIs.size() + inds.hfHadIs.size() == elements.size());
1280
1281
1282
1283
1284
1285 if (trackInBlock) {
1286
1287 std::multimap<double, unsigned> sortedTracks;
1288 std::multimap<double, unsigned> sortedTracksActive;
1289
1290 std::multimap<unsigned, std::pair<double, unsigned>> associatedHfEms;
1291
1292 std::multimap<double, std::pair<unsigned, double>> hfemSatellites;
1293
1294
1295
1296
1297 for (unsigned iHfHad : inds.hfHadIs) {
1298 PFBlockElement::Type type = elements[iHfHad].type();
1299 assert(type == PFBlockElement::HFHAD);
1300
1301 PFClusterRef hclusterRef = elements[iHfHad].clusterRef();
1302 assert(!hclusterRef.isNull());
1303
1304 sortedTracks.clear();
1305 sortedTracksActive.clear();
1306 associatedHfEms.clear();
1307 hfemSatellites.clear();
1308
1309
1310 block.associatedElements(
1311 iHfHad, linkData, sortedTracks, reco::PFBlockElement::TRACK, reco::PFBlock::LINKTEST_ALL);
1312
1313 LogTrace("PFAlgo|createCandidatesHF") << "elements[" << iHfHad << "]=" << elements[iHfHad];
1314
1315 if (sortedTracks.empty()) {
1316 LogTrace("PFAlgo|createCandidatesHCF") << "\tno associated tracks, keep for later";
1317 continue;
1318 }
1319
1320
1321 active[iHfHad] = false;
1322
1323 LogTrace("PFAlgo|createCandidatesHF") << sortedTracks.size() << " associated tracks:";
1324
1325 double totalChargedMomentum = 0.;
1326 double sumpError2 = 0.;
1327
1328
1329
1330
1331 for (auto const& trk : sortedTracks) {
1332 unsigned iTrack = trk.second;
1333
1334
1335 if (!active[iTrack])
1336 continue;
1337
1338 PFBlockElement::Type type = elements[iTrack].type();
1339 assert(type == reco::PFBlockElement::TRACK);
1340
1341 reco::TrackRef trackRef = elements[iTrack].trackRef();
1342 assert(!trackRef.isNull());
1343
1344
1345 double trackMomentum = trackRef->p();
1346 totalChargedMomentum += trackMomentum;
1347
1348
1349 double dp = trackRef->qoverpError() * trackMomentum * trackMomentum;
1350 sumpError2 += dp * dp;
1351
1352
1353 sortedTracksActive.emplace(trk);
1354
1355
1356 std::multimap<double, unsigned> sortedHfEms;
1357 block.associatedElements(
1358 iTrack, linkData, sortedHfEms, reco::PFBlockElement::HFEM, reco::PFBlock::LINKTEST_ALL);
1359
1360 LogTrace("PFAlgo|createCandidatesHF") << "number of HfEm elements linked to this track: " << sortedHfEms.size();
1361
1362 bool connectedToHfEm = false;
1363
1364
1365
1366
1367 for (auto const& hfem : sortedHfEms) {
1368 unsigned iHfEm = hfem.second;
1369 double dist = hfem.first;
1370
1371
1372 if (!active[iHfEm]) {
1373 LogTrace("PFAlgo|createCandidatesHF") << "cluster locked";
1374 continue;
1375 }
1376
1377
1378 PFBlockElement::Type type = elements[iHfEm].type();
1379 assert(type == PFBlockElement::HFEM);
1380 PFClusterRef eclusterRef = elements[iHfEm].clusterRef();
1381 assert(!eclusterRef.isNull());
1382
1383
1384 std::multimap<double, unsigned> sortedTracksHfEm;
1385 block.associatedElements(
1386 iHfEm, linkData, sortedTracksHfEm, reco::PFBlockElement::TRACK, reco::PFBlock::LINKTEST_ALL);
1387 unsigned jTrack = sortedTracksHfEm.begin()->second;
1388 if (jTrack != iTrack)
1389 continue;
1390
1391 double distHfEm = block.dist(jTrack, iHfEm, linkData, reco::PFBlock::LINKTEST_ALL);
1392 double hfemEnergy = eclusterRef->energy();
1393
1394 if (!connectedToHfEm) {
1395
1396 LogTrace("PFAlgo|createCandidatesHF") << "closest: " << elements[iHfEm];
1397 connectedToHfEm = true;
1398 std::pair<unsigned, double> satellite(iHfEm, hfemEnergy);
1399 hfemSatellites.emplace(-1., satellite);
1400
1401 } else {
1402
1403
1404 std::pair<unsigned, double> satellite(iHfEm, hfemEnergy);
1405 hfemSatellites.emplace(dist, satellite);
1406 }
1407
1408 std::pair<double, unsigned> associatedHfEm(distHfEm, iHfEm);
1409 associatedHfEms.emplace(iTrack, associatedHfEm);
1410
1411 }
1412 }
1413
1414
1415 double uncalibratedenergyHfHad = hclusterRef->energy();
1416 double energyHfHad = uncalibratedenergyHfHad;
1417 if (thepfEnergyCalibrationHF_.getcalibHF_use()) {
1418 energyHfHad = thepfEnergyCalibrationHF_.energyHad(
1419 uncalibratedenergyHfHad,
1420 hclusterRef->positionREP().Eta(),
1421 hclusterRef->positionREP().Phi());
1422 }
1423 double calibFactorHfHad = (uncalibratedenergyHfHad > 0.) ? energyHfHad / uncalibratedenergyHfHad : 1.;
1424
1425
1426 double energyHfEmTmp = 0.;
1427 double uncalibratedenergyHfEmTmp = 0.;
1428 double energyHfEm = 0.;
1429 double uncalibratedenergyHfEm = 0.;
1430
1431
1432 double caloResolution = hfEnergyResolution(totalChargedMomentum);
1433 caloResolution *= totalChargedMomentum;
1434 double totalError = sqrt(caloResolution * caloResolution + sumpError2);
1435 double nsigmaHFEM = nSigmaHFEM(totalChargedMomentum);
1436 double nsigmaHFHAD = nSigmaHFHAD(totalChargedMomentum);
1437
1438
1439 if (sortedTracksActive.empty()) {
1440
1441 std::multimap<double, unsigned> sortedHfEms;
1442 std::multimap<double, unsigned> sortedHfEmsActive;
1443 block.associatedElements(
1444 iHfHad, linkData, sortedHfEms, reco::PFBlockElement::HFEM, reco::PFBlock::LINKTEST_ALL);
1445
1446
1447
1448 if (!sortedHfEms.empty()) {
1449 for (auto const& hfem : sortedHfEms) {
1450 unsigned iHfEm = hfem.second;
1451
1452 if (!active[iHfEm])
1453 continue;
1454 sortedHfEmsActive.emplace(hfem);
1455 PFClusterRef eclusterRef = elements[iHfEm].clusterRef();
1456 assert(!eclusterRef.isNull());
1457 double hfemEnergy = eclusterRef->energy();
1458 uncalibratedenergyHfEm += hfemEnergy;
1459 energyHfEm = uncalibratedenergyHfEm;
1460 if (thepfEnergyCalibrationHF_.getcalibHF_use()) {
1461 energyHfEm = thepfEnergyCalibrationHF_.energyEmHad(
1462 uncalibratedenergyHfEm, 0.0, eclusterRef->positionREP().Eta(), eclusterRef->positionREP().Phi());
1463 energyHfHad = thepfEnergyCalibrationHF_.energyEmHad(
1464 0.0, uncalibratedenergyHfHad, hclusterRef->positionREP().Eta(), hclusterRef->positionREP().Phi());
1465 }
1466 }
1467 }
1468
1469
1470 unsigned tmpi = reconstructCluster(*hclusterRef, energyHfEm + energyHfHad);
1471 (*pfCandidates_)[tmpi].setHcalEnergy(uncalibratedenergyHfHad, energyHfHad);
1472 (*pfCandidates_)[tmpi].setEcalEnergy(uncalibratedenergyHfEm, energyHfEm);
1473 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iHfHad);
1474 for (auto const& hfem : sortedHfEmsActive) {
1475 unsigned iHfEm = hfem.second;
1476 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iHfEm);
1477 active[iHfEm] = false;
1478 }
1479
1480 }
1481
1482
1483
1484 else if ((energyHfHad - totalChargedMomentum) > nsigmaHFHAD * totalError) {
1485 assert(energyHfEm == 0.);
1486
1487 double energyHfHadExcess = max(energyHfHad - totalChargedMomentum, 0.);
1488 double uncalibratedenergyHfHadExcess = energyHfHadExcess / calibFactorHfHad;
1489 unsigned tmpi = reconstructCluster(*hclusterRef, energyHfHadExcess);
1490 (*pfCandidates_)[tmpi].setHcalEnergy(uncalibratedenergyHfHadExcess, energyHfHadExcess);
1491 (*pfCandidates_)[tmpi].setEcalEnergy(0., 0.);
1492 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iHfHad);
1493 energyHfHad = max(energyHfHad - energyHfHadExcess, 0.);
1494 uncalibratedenergyHfHad = max(uncalibratedenergyHfHad - uncalibratedenergyHfHadExcess, 0.);
1495 }
1496
1497
1498
1499
1500
1501 else {
1502 for (auto const& hfemSatellite : hfemSatellites) {
1503
1504 uncalibratedenergyHfEmTmp += std::get<1>(hfemSatellite.second);
1505 energyHfEmTmp = uncalibratedenergyHfEmTmp;
1506 double energyHfHadTmp = uncalibratedenergyHfHad;
1507 unsigned iHfEm = std::get<0>(hfemSatellite.second);
1508 PFClusterRef eclusterRef = elements[iHfEm].clusterRef();
1509 assert(!eclusterRef.isNull());
1510 if (thepfEnergyCalibrationHF_.getcalibHF_use()) {
1511 energyHfEmTmp = thepfEnergyCalibrationHF_.energyEmHad(
1512 uncalibratedenergyHfEmTmp, 0.0, eclusterRef->positionREP().Eta(), eclusterRef->positionREP().Phi());
1513 energyHfHadTmp = thepfEnergyCalibrationHF_.energyEmHad(
1514 0.0, uncalibratedenergyHfHad, hclusterRef->positionREP().Eta(), hclusterRef->positionREP().Phi());
1515 }
1516
1517 double caloEnergyTmp = energyHfEmTmp + energyHfHadTmp;
1518 double calibFactorHfEm = (uncalibratedenergyHfEmTmp > 0.) ? energyHfEmTmp / uncalibratedenergyHfEmTmp : 1.;
1519
1520
1521
1522 if (hfemSatellite.first < 0. || caloEnergyTmp < totalChargedMomentum) {
1523 LogTrace("PFAlgo|createCandidatesHF")
1524 << "\t\t\tactive, adding " << std::get<1>(hfemSatellite.second) << " to HFEM energy, and locking";
1525 active[std::get<0>(hfemSatellite.second)] = false;
1526
1527 if (hfemSatellite.first < 0. && (caloEnergyTmp - totalChargedMomentum) > nsigmaHFEM * totalError) {
1528
1529 double energyHfEmExcess = max(caloEnergyTmp - totalChargedMomentum, 0.);
1530 double uncalibratedenergyHfEmExcess = energyHfEmExcess / calibFactorHfEm;
1531 unsigned tmpi = reconstructCluster(*eclusterRef, energyHfEmExcess);
1532 (*pfCandidates_)[tmpi].setEcalEnergy(uncalibratedenergyHfEmExcess, energyHfEmExcess);
1533 (*pfCandidates_)[tmpi].setHcalEnergy(0, 0.);
1534 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iHfEm);
1535 energyHfEmTmp = max(energyHfEmTmp - energyHfEmExcess, 0.);
1536 uncalibratedenergyHfEmTmp = max(uncalibratedenergyHfEmTmp - uncalibratedenergyHfEmExcess, 0.);
1537 }
1538 energyHfEm = energyHfEmTmp;
1539 uncalibratedenergyHfEm = uncalibratedenergyHfEmTmp;
1540 energyHfHad = energyHfHadTmp;
1541 continue;
1542 }
1543 break;
1544 }
1545 }
1546
1547
1548
1549
1550 for (auto const& trk : sortedTracksActive) {
1551 unsigned iTrack = trk.second;
1552
1553
1554 reco::TrackRef trackRef = elements[iTrack].trackRef();
1555 assert(!trackRef.isNull());
1556
1557
1558
1559
1560 unsigned tmpi = reconstructTrack(elements[iTrack]);
1561 active[iTrack] = false;
1562 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iHfHad);
1563 auto myHfEms = associatedHfEms.equal_range(iTrack);
1564 for (auto ii = myHfEms.first; ii != myHfEms.second; ++ii) {
1565 unsigned iHfEm = ii->second.second;
1566 if (active[iHfEm])
1567 continue;
1568 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iHfEm);
1569 }
1570 double frac = 0.;
1571 if (totalChargedMomentum)
1572 frac = trackRef->p() / totalChargedMomentum;
1573 (*pfCandidates_)[tmpi].setEcalEnergy(uncalibratedenergyHfEm * frac, energyHfEm * frac);
1574 (*pfCandidates_)[tmpi].setHcalEnergy(uncalibratedenergyHfHad * frac, energyHfHad * frac);
1575
1576 }
1577
1578 }
1579
1580
1581
1582
1583 for (unsigned iHfEm = 0; iHfEm < elements.size(); iHfEm++) {
1584 PFBlockElement::Type type = elements[iHfEm].type();
1585 if (type == PFBlockElement::HFEM && active[iHfEm]) {
1586 reco::PFClusterRef eclusterRef = elements[iHfEm].clusterRef();
1587 double energyHF = 0.;
1588 double uncalibratedenergyHF = 0.;
1589 unsigned tmpi = 0;
1590
1591 energyHF = eclusterRef->energy();
1592 uncalibratedenergyHF = energyHF;
1593 if (thepfEnergyCalibrationHF_.getcalibHF_use()) {
1594 energyHF = thepfEnergyCalibrationHF_.energyEm(
1595 uncalibratedenergyHF, eclusterRef->positionREP().Eta(), eclusterRef->positionREP().Phi());
1596 }
1597 tmpi = reconstructCluster(*eclusterRef, energyHF);
1598 (*pfCandidates_)[tmpi].setEcalEnergy(uncalibratedenergyHF, energyHF);
1599 (*pfCandidates_)[tmpi].setHcalEnergy(0., 0.);
1600 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iHfEm);
1601 active[iHfEm] = false;
1602 LogTrace("PFAlgo|createCandidatesHF") << "HF EM alone from blocks with tracks! " << energyHF;
1603 }
1604 }
1605
1606 }
1607
1608
1609
1610
1611
1612 else if (elements.size() == 1) {
1613
1614 reco::PFClusterRef clusterRef = elements[0].clusterRef();
1615 double energyHF = 0.;
1616 double uncalibratedenergyHF = 0.;
1617 unsigned tmpi = 0;
1618 switch (clusterRef->layer()) {
1619 case PFLayer::HF_EM:
1620
1621 energyHF = clusterRef->energy();
1622 uncalibratedenergyHF = energyHF;
1623 if (thepfEnergyCalibrationHF_.getcalibHF_use()) {
1624 energyHF = thepfEnergyCalibrationHF_.energyEm(
1625 uncalibratedenergyHF, clusterRef->positionREP().Eta(), clusterRef->positionREP().Phi());
1626 }
1627 tmpi = reconstructCluster(*clusterRef, energyHF);
1628 (*pfCandidates_)[tmpi].setEcalEnergy(uncalibratedenergyHF, energyHF);
1629 (*pfCandidates_)[tmpi].setHcalEnergy(0., 0.);
1630 (*pfCandidates_)[tmpi].setHoEnergy(0., 0.);
1631 (*pfCandidates_)[tmpi].setPs1Energy(0.);
1632 (*pfCandidates_)[tmpi].setPs2Energy(0.);
1633 (*pfCandidates_)[tmpi].addElementInBlock(blockref, inds.hfEmIs[0]);
1634 LogTrace("PFAlgo|createCandidatesHF") << "HF EM alone ! " << energyHF;
1635 break;
1636 case PFLayer::HF_HAD:
1637
1638 energyHF = clusterRef->energy();
1639 uncalibratedenergyHF = energyHF;
1640 if (thepfEnergyCalibrationHF_.getcalibHF_use()) {
1641 energyHF = thepfEnergyCalibrationHF_.energyHad(
1642 uncalibratedenergyHF, clusterRef->positionREP().Eta(), clusterRef->positionREP().Phi());
1643 }
1644 tmpi = reconstructCluster(*clusterRef, energyHF);
1645 (*pfCandidates_)[tmpi].setHcalEnergy(uncalibratedenergyHF, energyHF);
1646 (*pfCandidates_)[tmpi].setEcalEnergy(0., 0.);
1647 (*pfCandidates_)[tmpi].setHoEnergy(0., 0.);
1648 (*pfCandidates_)[tmpi].setPs1Energy(0.);
1649 (*pfCandidates_)[tmpi].setPs2Energy(0.);
1650 (*pfCandidates_)[tmpi].addElementInBlock(blockref, inds.hfHadIs[0]);
1651 LogTrace("PFAlgo|createCandidatesHF") << "HF Had alone ! " << energyHF;
1652 break;
1653 default:
1654 assert(0);
1655 }
1656 } else if (elements.size() == 2) {
1657
1658 reco::PFClusterRef c0 = elements[0].clusterRef();
1659 reco::PFClusterRef c1 = elements[1].clusterRef();
1660
1661 reco::PFClusterRef cem = (c0->layer() == PFLayer::HF_EM ? c0 : c1);
1662 reco::PFClusterRef chad = (c1->layer() == PFLayer::HF_HAD ? c1 : c0);
1663
1664 if (cem->layer() != PFLayer::HF_EM || chad->layer() != PFLayer::HF_HAD) {
1665 edm::LogError("PFAlgo::createCandidatesHF") << "Error: 2 elements, but not 1 HFEM and 1 HFHAD";
1666 edm::LogError("PFAlgo::createCandidatesHF") << block;
1667 assert(0);
1668
1669
1670 }
1671
1672 double energyHfEm = cem->energy();
1673 double energyHfHad = chad->energy();
1674 double uncalibratedenergyHfEm = energyHfEm;
1675 double uncalibratedenergyHfHad = energyHfHad;
1676 if (thepfEnergyCalibrationHF_.getcalibHF_use()) {
1677 energyHfEm = thepfEnergyCalibrationHF_.energyEmHad(
1678 uncalibratedenergyHfEm, 0.0, c0->positionREP().Eta(), c0->positionREP().Phi());
1679 energyHfHad = thepfEnergyCalibrationHF_.energyEmHad(
1680 0.0, uncalibratedenergyHfHad, c1->positionREP().Eta(), c1->positionREP().Phi());
1681 }
1682 auto& cand = (*pfCandidates_)[reconstructCluster(*chad, energyHfEm + energyHfHad)];
1683 cand.setEcalEnergy(uncalibratedenergyHfEm, energyHfEm);
1684 cand.setHcalEnergy(uncalibratedenergyHfHad, energyHfHad);
1685 cand.setHoEnergy(0., 0.);
1686 cand.setPs1Energy(0.);
1687 cand.setPs2Energy(0.);
1688 cand.addElementInBlock(blockref, inds.hfEmIs[0]);
1689 cand.addElementInBlock(blockref, inds.hfHadIs[0]);
1690 LogTrace("PFAlgo|createCandidatesHF") << "HF EM+HAD found ! " << energyHfEm << " " << energyHfHad;
1691 } else {
1692
1693 edm::LogWarning("PFAlgo::createCandidatesHF")
1694 << "Warning: HF, but n elem different from 1 or 2 or >=3 or !trackIs.empty()";
1695 edm::LogWarning("PFAlgo::createCandidatesHF") << block;
1696 }
1697 LogTrace("PFAlgo|createCandidateHF") << "end of function PFAlgo::createCandidateHF";
1698 }
1699
1700 void PFAlgo::createCandidatesHCAL(const reco::PFBlock& block,
1701 reco::PFBlock::LinkData& linkData,
1702 const edm::OwnVector<reco::PFBlockElement>& elements,
1703 std::vector<bool>& active,
1704 const reco::PFBlockRef& blockref,
1705 ElementIndices& inds,
1706 std::vector<bool>& deadArea) {
1707 LogTrace("PFAlgo|createCandidatesHCAL")
1708 << "start of function PFAlgo::createCandidatesHCAL, inds.hcalIs.size()=" << inds.hcalIs.size();
1709
1710
1711
1712 for (unsigned iHcal : inds.hcalIs) {
1713 PFBlockElement::Type type = elements[iHcal].type();
1714
1715 assert(type == PFBlockElement::HCAL);
1716
1717 LogTrace("PFAlgo|createCandidatesHCAL") << "elements[" << iHcal << "]=" << elements[iHcal];
1718
1719
1720 std::multimap<double, unsigned> sortedTracks;
1721 block.associatedElements(iHcal, linkData, sortedTracks, reco::PFBlockElement::TRACK, reco::PFBlock::LINKTEST_ALL);
1722
1723 std::multimap<unsigned, std::pair<double, unsigned>> associatedEcals;
1724
1725 std::map<unsigned, std::pair<double, double>> associatedPSs;
1726
1727 std::multimap<double, std::pair<unsigned, bool>> associatedTracks;
1728
1729
1730 std::multimap<double, std::tuple<unsigned, ::math::XYZVector, double>>
1731 ecalSatellites;
1732 std::tuple<unsigned, ::math::XYZVector, double> fakeSatellite(iHcal, ::math::XYZVector(0., 0., 0.), 1.);
1733 ecalSatellites.emplace(-1., fakeSatellite);
1734
1735 std::multimap<unsigned, std::pair<double, unsigned>> associatedHOs;
1736
1737 PFClusterRef hclusterref = elements[iHcal].clusterRef();
1738 assert(!hclusterref.isNull());
1739
1740
1741
1742
1743
1744
1745
1746
1747 if (sortedTracks.empty()) {
1748 LogTrace("PFAlgo|createCandidatesHCAL") << "\tno associated tracks, keep for later";
1749 continue;
1750 }
1751
1752
1753 active[iHcal] = false;
1754
1755
1756
1757
1758
1759
1760 LogTrace("PFAlgo|createCandidatesHCAL") << sortedTracks.size() << " associated tracks:";
1761
1762 double totalChargedMomentum = 0;
1763 double sumpError2 = 0.;
1764 double totalHO = 0.;
1765 double totalEcal = 0.;
1766 double totalEcalEGMCalib = 0.;
1767 double totalHcal = hclusterref->energy();
1768 vector<double> hcalP;
1769 vector<double> hcalDP;
1770 vector<unsigned> tkIs;
1771 double maxDPovP = -9999.;
1772
1773
1774 vector<unsigned> chargedHadronsIndices;
1775 vector<unsigned> chargedHadronsInBlock;
1776 double mergedNeutralHadronEnergy = 0;
1777 double mergedPhotonEnergy = 0;
1778 double muonHCALEnergy = 0.;
1779 double muonECALEnergy = 0.;
1780 double muonHCALError = 0.;
1781 double muonECALError = 0.;
1782 unsigned nMuons = 0;
1783
1784 ::math::XYZVector photonAtECAL(0., 0., 0.);
1785 std::vector<std::tuple<unsigned, ::math::XYZVector, double>>
1786 ecalClusters;
1787 double sumEcalClusters = 0;
1788 ::math::XYZVector hadronDirection(
1789 hclusterref->position().X(), hclusterref->position().Y(), hclusterref->position().Z());
1790 hadronDirection = hadronDirection.Unit();
1791 ::math::XYZVector hadronAtECAL = totalHcal * hadronDirection;
1792
1793
1794 for (auto const& trk : sortedTracks) {
1795 unsigned iTrack = trk.second;
1796
1797
1798 if (!active[iTrack])
1799 continue;
1800
1801 PFBlockElement::Type type = elements[iTrack].type();
1802 assert(type == reco::PFBlockElement::TRACK);
1803
1804 reco::TrackRef trackRef = elements[iTrack].trackRef();
1805 assert(!trackRef.isNull());
1806
1807
1808 const ::math::XYZPointF& chargedPosition =
1809 dynamic_cast<const reco::PFBlockElementTrack*>(&elements[iTrack])->positionAtECALEntrance();
1810 ::math::XYZVector chargedDirection(chargedPosition.X(), chargedPosition.Y(), chargedPosition.Z());
1811 chargedDirection = chargedDirection.Unit();
1812
1813
1814 std::multimap<double, unsigned> sortedEcals;
1815 block.associatedElements(iTrack, linkData, sortedEcals, reco::PFBlockElement::ECAL, reco::PFBlock::LINKTEST_ALL);
1816
1817 LogTrace("PFAlgo|createCandidatesHCAL") << "number of Ecal elements linked to this track: " << sortedEcals.size();
1818
1819
1820 std::multimap<double, unsigned> sortedHOs;
1821 if (useHO_) {
1822 block.associatedElements(iTrack, linkData, sortedHOs, reco::PFBlockElement::HO, reco::PFBlock::LINKTEST_ALL);
1823 }
1824 LogTrace("PFAlgo|createCandidatesHCAL")
1825 << "PFAlgo : number of HO elements linked to this track: " << sortedHOs.size();
1826
1827
1828 reco::MuonRef muonRef = elements[iTrack].muonRef();
1829
1830 bool thisIsAMuon = PFMuonAlgo::isMuon(elements[iTrack]);
1831 bool thisIsAnIsolatedMuon = PFMuonAlgo::isIsolatedMuon(elements[iTrack]);
1832 bool thisIsALooseMuon = false;
1833
1834 if (!thisIsAMuon) {
1835 thisIsALooseMuon = PFMuonAlgo::isLooseMuon(elements[iTrack]);
1836 }
1837
1838 if (thisIsAMuon) {
1839 LogTrace("PFAlgo|createCandidatesHCAL") << "This track is identified as a muon - remove it from the stack";
1840 LogTrace("PFAlgo|createCandidatesHCAL") << elements[iTrack];
1841
1842
1843
1844 unsigned tmpi = reconstructTrack(elements[iTrack]);
1845
1846 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iTrack);
1847 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iHcal);
1848 double muonHcal = std::min(muonHCAL_[0] + muonHCAL_[1], totalHcal);
1849
1850
1851
1852
1853 bool letMuonEatCaloEnergy = false;
1854
1855 if (thisIsAnIsolatedMuon) {
1856
1857 double totalCaloEnergy = totalHcal / 1.30;
1858 unsigned iEcal = 0;
1859 if (!sortedEcals.empty()) {
1860 iEcal = sortedEcals.begin()->second;
1861 PFClusterRef eclusterref = elements[iEcal].clusterRef();
1862 totalCaloEnergy += eclusterref->energy();
1863 }
1864
1865 if (useHO_) {
1866
1867 unsigned iHO = 0;
1868 if (!sortedHOs.empty()) {
1869 iHO = sortedHOs.begin()->second;
1870 PFClusterRef eclusterref = elements[iHO].clusterRef();
1871 totalCaloEnergy += eclusterref->energy() / 1.30;
1872 }
1873 }
1874
1875 if ((pfCandidates_->back()).p() > totalCaloEnergy)
1876 letMuonEatCaloEnergy = true;
1877 }
1878
1879 if (letMuonEatCaloEnergy)
1880 muonHcal = totalHcal;
1881 double muonEcal = 0.;
1882 unsigned iEcal = 0;
1883 if (!sortedEcals.empty()) {
1884 iEcal = sortedEcals.begin()->second;
1885 PFClusterRef eclusterref = elements[iEcal].clusterRef();
1886 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iEcal);
1887 muonEcal = std::min(muonECAL_[0] + muonECAL_[1], eclusterref->energy());
1888 if (letMuonEatCaloEnergy)
1889 muonEcal = eclusterref->energy();
1890
1891 if (eclusterref->energy() - muonEcal < 0.2)
1892 active[iEcal] = false;
1893 (*pfCandidates_)[tmpi].setEcalEnergy(eclusterref->energy(), muonEcal);
1894 }
1895 unsigned iHO = 0;
1896 double muonHO = 0.;
1897 if (useHO_) {
1898 if (!sortedHOs.empty()) {
1899 iHO = sortedHOs.begin()->second;
1900 PFClusterRef hoclusterref = elements[iHO].clusterRef();
1901 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iHO);
1902 muonHO = std::min(muonHO_[0] + muonHO_[1], hoclusterref->energy());
1903 if (letMuonEatCaloEnergy)
1904 muonHO = hoclusterref->energy();
1905
1906 if (hoclusterref->energy() - muonHO < 0.2)
1907 active[iHO] = false;
1908 (*pfCandidates_)[tmpi].setHcalEnergy(totalHcal, muonHcal);
1909 (*pfCandidates_)[tmpi].setHoEnergy(hoclusterref->energy(), muonHO);
1910 }
1911 } else {
1912 (*pfCandidates_)[tmpi].setHcalEnergy(totalHcal, muonHcal);
1913 }
1914 setHcalDepthInfo((*pfCandidates_)[tmpi], *hclusterref);
1915
1916 if (letMuonEatCaloEnergy) {
1917 muonHCALEnergy += totalHcal;
1918 if (useHO_)
1919 muonHCALEnergy += muonHO;
1920 muonHCALError += 0.;
1921 muonECALEnergy += muonEcal;
1922 muonECALError += 0.;
1923 photonAtECAL -= muonEcal * chargedDirection;
1924 hadronAtECAL -= totalHcal * chargedDirection;
1925 if (!sortedEcals.empty())
1926 active[iEcal] = false;
1927 active[iHcal] = false;
1928 if (useHO_ && !sortedHOs.empty())
1929 active[iHO] = false;
1930 } else {
1931
1932 muonHCALEnergy += muonHCAL_[0];
1933 muonHCALError += muonHCAL_[1] * muonHCAL_[1];
1934 if (muonHO > 0.) {
1935 muonHCALEnergy += muonHO_[0];
1936 muonHCALError += muonHO_[1] * muonHO_[1];
1937 }
1938 muonECALEnergy += muonECAL_[0];
1939 muonECALError += muonECAL_[1] * muonECAL_[1];
1940
1941 photonAtECAL -= muonECAL_[0] * chargedDirection;
1942 hadronAtECAL -= muonHCAL_[0] * chargedDirection;
1943 }
1944
1945
1946 active[iTrack] = false;
1947
1948 continue;
1949 }
1950
1951
1952
1953 LogTrace("PFAlgo|createCandidatesHCAL") << "elements[iTrack=" << iTrack << "]=" << elements[iTrack];
1954
1955
1956 double trackMomentum = trackRef->p();
1957 totalChargedMomentum += trackMomentum;
1958
1959
1960
1961 if (thisIsALooseMuon && !thisIsAMuon)
1962 nMuons += 1;
1963
1964
1965
1966
1967 double dpt = trackRef->ptError();
1968 double blowError = PFTrackAlgoTools::errorScale(trackRef->algo(), factors45_);
1969
1970 bool isPrimaryOrSecondary = isFromSecInt(elements[iTrack], "all");
1971
1972 if (isPrimaryOrSecondary)
1973 blowError = 1.;
1974
1975 std::pair<unsigned, bool> tkmuon(iTrack, thisIsALooseMuon);
1976 associatedTracks.emplace(-dpt * blowError, tkmuon);
1977
1978
1979 double dp = trackRef->qoverpError() * trackMomentum * trackMomentum;
1980 sumpError2 += dp * dp;
1981
1982 bool connectedToEcal = false;
1983 if (!sortedEcals.empty()) {
1984
1985
1986 for (auto const& ecal : sortedEcals) {
1987 unsigned iEcal = ecal.second;
1988 double dist = ecal.first;
1989
1990
1991 if (!active[iEcal]) {
1992 LogTrace("PFAlgo|createCandidatesHCAL") << "cluster locked";
1993 continue;
1994 }
1995
1996
1997 PFBlockElement::Type type = elements[iEcal].type();
1998 assert(type == PFBlockElement::ECAL);
1999 PFClusterRef eclusterref = elements[iEcal].clusterRef();
2000 assert(!eclusterref.isNull());
2001
2002
2003 std::multimap<double, unsigned> sortedTracksEcal;
2004 block.associatedElements(
2005 iEcal, linkData, sortedTracksEcal, reco::PFBlockElement::TRACK, reco::PFBlock::LINKTEST_ALL);
2006 unsigned jTrack = sortedTracksEcal.begin()->second;
2007 if (jTrack != iTrack)
2008 continue;
2009
2010 double distEcal = block.dist(jTrack, iEcal, linkData, reco::PFBlock::LINKTEST_ALL);
2011
2012 float ecalEnergyCalibrated = eclusterref->correctedEnergy();
2013 float ecalEnergy = eclusterref->energy();
2014 ::math::XYZVector photonDirection(
2015 eclusterref->position().X(), eclusterref->position().Y(), eclusterref->position().Z());
2016 photonDirection = photonDirection.Unit();
2017
2018 if (!connectedToEcal) {
2019
2020 LogTrace("PFAlgo|createCandidatesHCAL") << "closest: " << elements[iEcal];
2021
2022 connectedToEcal = true;
2023
2024
2025
2026
2027 double ecalCalibFactor = (ecalEnergy > 1E-9) ? ecalEnergyCalibrated / ecalEnergy : 1.;
2028 std::tuple<unsigned, ::math::XYZVector, double> satellite(
2029 iEcal, ecalEnergy * photonDirection, ecalCalibFactor);
2030 ecalSatellites.emplace(-1., satellite);
2031
2032 } else {
2033
2034
2035 double ecalCalibFactor = (ecalEnergy > 1E-9) ? ecalEnergyCalibrated / ecalEnergy : 1.;
2036 std::tuple<unsigned, ::math::XYZVector, double> satellite(
2037 iEcal, ecalEnergy * photonDirection, ecalCalibFactor);
2038 ecalSatellites.emplace(dist, satellite);
2039 }
2040
2041 std::pair<double, unsigned> associatedEcal(distEcal, iEcal);
2042 associatedEcals.emplace(iTrack, associatedEcal);
2043
2044 }
2045 }
2046
2047 if (useHO_ && !sortedHOs.empty()) {
2048
2049
2050 for (auto const& ho : sortedHOs) {
2051 unsigned iHO = ho.second;
2052 double distHO = ho.first;
2053
2054
2055 if (!active[iHO]) {
2056 LogTrace("PFAlgo|createCandidatesHCAL") << "cluster locked";
2057 continue;
2058 }
2059
2060
2061 PFBlockElement::Type type = elements[iHO].type();
2062 assert(type == PFBlockElement::HO);
2063 PFClusterRef hoclusterref = elements[iHO].clusterRef();
2064 assert(!hoclusterref.isNull());
2065
2066
2067 std::multimap<double, unsigned> sortedTracksHO;
2068 block.associatedElements(
2069 iHO, linkData, sortedTracksHO, reco::PFBlockElement::TRACK, reco::PFBlock::LINKTEST_ALL);
2070 unsigned jTrack = sortedTracksHO.begin()->second;
2071 if (jTrack != iTrack)
2072 continue;
2073
2074
2075
2076
2077
2078
2079
2080 totalHO += hoclusterref->energy();
2081 active[iHO] = false;
2082
2083 std::pair<double, unsigned> associatedHO(distHO, iHO);
2084 associatedHOs.emplace(iTrack, associatedHO);
2085
2086 }
2087 }
2088
2089 }
2090
2091
2092 totalHcal += totalHO;
2093
2094
2095
2096 double caloEnergy = 0.;
2097 double slopeEcal = 1.0;
2098 double calibEcal = 0.;
2099 double calibHcal = 0.;
2100 hadronDirection = hadronAtECAL.Unit();
2101
2102
2103 double caloResolution = neutralHadronEnergyResolution(totalChargedMomentum, hclusterref->positionREP().Eta());
2104 caloResolution *= totalChargedMomentum;
2105
2106 caloResolution = std::sqrt(caloResolution * caloResolution + muonHCALError + muonECALError);
2107 totalEcal -= std::min(totalEcal, muonECALEnergy);
2108 totalEcalEGMCalib -= std::min(totalEcalEGMCalib, muonECALEnergy);
2109 totalHcal -= std::min(totalHcal, muonHCALEnergy);
2110 if (totalEcal < 1E-9)
2111 photonAtECAL = ::math::XYZVector(0., 0., 0.);
2112 if (totalHcal < 1E-9)
2113 hadronAtECAL = ::math::XYZVector(0., 0., 0.);
2114
2115
2116
2117
2118
2119 for (auto const& ecalSatellite : ecalSatellites) {
2120
2121 double previousCalibEcal = calibEcal;
2122 double previousCalibHcal = calibHcal;
2123 double previousCaloEnergy = caloEnergy;
2124 double previousSlopeEcal = slopeEcal;
2125 ::math::XYZVector previousHadronAtECAL = hadronAtECAL;
2126
2127 totalEcal +=
2128 sqrt(std::get<1>(ecalSatellite.second).Mag2());
2129 totalEcalEGMCalib += sqrt(std::get<1>(ecalSatellite.second).Mag2()) *
2130 std::get<2>(ecalSatellite.second);
2131 photonAtECAL += std::get<1>(ecalSatellite.second) *
2132 std::get<2>(ecalSatellite.second);
2133 calibEcal = std::max(0., totalEcal);
2134 calibHcal = std::max(0., totalHcal);
2135 hadronAtECAL = calibHcal * hadronDirection;
2136
2137 calibration_.energyEmHad(totalChargedMomentum,
2138 calibEcal,
2139 calibHcal,
2140 hclusterref->positionREP().Eta(),
2141 hclusterref->positionREP().Phi());
2142 caloEnergy = calibEcal + calibHcal;
2143 if (totalEcal > 0.)
2144 slopeEcal = calibEcal / totalEcal;
2145
2146 hadronAtECAL = calibHcal * hadronDirection;
2147
2148
2149
2150 if (ecalSatellite.first < 0. || caloEnergy - totalChargedMomentum <= 0.) {
2151 LogTrace("PFAlgo|createCandidatesHCAL")
2152 << "\t\t\tactive, adding " << std::get<1>(ecalSatellite.second) << " to ECAL energy, and locking";
2153 active[std::get<0>(ecalSatellite.second)] = false;
2154 double clusterEnergy =
2155 sqrt(std::get<1>(ecalSatellite.second).Mag2()) *
2156 std::get<2>(ecalSatellite.second);
2157 if (clusterEnergy > 50) {
2158 ecalClusters.push_back(ecalSatellite.second);
2159 sumEcalClusters += clusterEnergy;
2160 }
2161 continue;
2162 }
2163
2164
2165
2166 totalEcal -= sqrt(std::get<1>(ecalSatellite.second).Mag2());
2167 totalEcalEGMCalib -= sqrt(std::get<1>(ecalSatellite.second).Mag2()) * std::get<2>(ecalSatellite.second);
2168 photonAtECAL -= std::get<1>(ecalSatellite.second) * std::get<2>(ecalSatellite.second);
2169 calibEcal = previousCalibEcal;
2170 calibHcal = previousCalibHcal;
2171 hadronAtECAL = previousHadronAtECAL;
2172 slopeEcal = previousSlopeEcal;
2173 caloEnergy = previousCaloEnergy;
2174
2175 break;
2176 }
2177
2178
2179 assert(caloEnergy >= 0);
2180
2181
2182
2183
2184
2185
2186
2187
2188
2189
2190
2191
2192
2193 if (totalChargedMomentum - caloEnergy > nSigmaTRACK_ * caloResolution) {
2194
2195 if (nMuons > 0) {
2196 for (auto const& trk : associatedTracks) {
2197
2198 if (!trk.second.second)
2199 continue;
2200
2201 const unsigned int iTrack = trk.second.first;
2202
2203 if (!active[iTrack])
2204 continue;
2205
2206 const double trackMomentum = elements[trk.second.first].trackRef()->p();
2207
2208
2209 std::multimap<double, unsigned> sortedEcals;
2210 block.associatedElements(
2211 iTrack, linkData, sortedEcals, reco::PFBlockElement::ECAL, reco::PFBlock::LINKTEST_ALL);
2212 std::multimap<double, unsigned> sortedHOs;
2213 block.associatedElements(iTrack, linkData, sortedHOs, reco::PFBlockElement::HO, reco::PFBlock::LINKTEST_ALL);
2214
2215
2216 auto& muon = (*pfCandidates_)[reconstructTrack(elements[iTrack], true)];
2217
2218 muon.addElementInBlock(blockref, iTrack);
2219 muon.addElementInBlock(blockref, iHcal);
2220 const double muonHcal = std::min(muonHCAL_[0] + muonHCAL_[1], totalHcal - totalHO);
2221 double muonHO = 0.;
2222 muon.setHcalEnergy(totalHcal, muonHcal);
2223 if (!sortedEcals.empty()) {
2224 const unsigned int iEcal = sortedEcals.begin()->second;
2225 PFClusterRef eclusterref = elements[iEcal].clusterRef();
2226 muon.addElementInBlock(blockref, iEcal);
2227 const double muonEcal = std::min(muonECAL_[0] + muonECAL_[1], eclusterref->energy());
2228 muon.setEcalEnergy(eclusterref->energy(), muonEcal);
2229 }
2230 if (useHO_ && !sortedHOs.empty()) {
2231 const unsigned int iHO = sortedHOs.begin()->second;
2232 PFClusterRef hoclusterref = elements[iHO].clusterRef();
2233 muon.addElementInBlock(blockref, iHO);
2234 muonHO = std::min(muonHO_[0] + muonHO_[1], hoclusterref->energy());
2235 muon.setHcalEnergy(max(totalHcal - totalHO, 0.0), muonHcal);
2236 muon.setHoEnergy(hoclusterref->energy(), muonHO);
2237 }
2238 setHcalDepthInfo(muon, *hclusterref);
2239
2240 const ::math::XYZPointF& chargedPosition =
2241 dynamic_cast<const reco::PFBlockElementTrack*>(&elements[trk.second.first])->positionAtECALEntrance();
2242 ::math::XYZVector chargedDirection(chargedPosition.X(), chargedPosition.Y(), chargedPosition.Z());
2243 chargedDirection = chargedDirection.Unit();
2244 totalChargedMomentum -= trackMomentum;
2245
2246 if (totalEcal > 0.)
2247 calibEcal -= std::min(calibEcal, muonECAL_[0] * calibEcal / totalEcal);
2248 if (totalHcal > 0.)
2249 calibHcal -= std::min(calibHcal, muonHCAL_[0] * calibHcal / totalHcal);
2250 totalEcal -= std::min(totalEcal, muonECAL_[0]);
2251 totalHcal -= std::min(totalHcal, muonHCAL_[0]);
2252 if (totalEcal > muonECAL_[0])
2253 photonAtECAL -= muonECAL_[0] * chargedDirection;
2254 if (totalHcal > muonHCAL_[0])
2255 hadronAtECAL -= muonHCAL_[0] * calibHcal / totalHcal * chargedDirection;
2256 caloEnergy = calibEcal + calibHcal;
2257 muonHCALEnergy += muonHCAL_[0];
2258 muonHCALError += muonHCAL_[1] * muonHCAL_[1];
2259 if (muonHO > 0.) {
2260 muonHCALEnergy += muonHO_[0];
2261 muonHCALError += muonHO_[1] * muonHO_[1];
2262 if (totalHcal > 0.) {
2263 calibHcal -= std::min(calibHcal, muonHO_[0] * calibHcal / totalHcal);
2264 totalHcal -= std::min(totalHcal, muonHO_[0]);
2265 }
2266 }
2267 muonECALEnergy += muonECAL_[0];
2268 muonECALError += muonECAL_[1] * muonECAL_[1];
2269 active[iTrack] = false;
2270
2271
2272
2273 }
2274
2275 caloResolution = neutralHadronEnergyResolution(totalChargedMomentum, hclusterref->positionREP().Eta());
2276 caloResolution *= totalChargedMomentum;
2277 caloResolution = std::sqrt(caloResolution * caloResolution + muonHCALError + muonECALError);
2278 }
2279 }
2280
2281 #ifdef EDM_ML_DEBUG
2282 LogTrace("PFAlgo|createCandidatesHCAL") << "\tBefore Cleaning ";
2283 LogTrace("PFAlgo|createCandidatesHCAL") << "\tCompare Calo Energy to total charged momentum ";
2284 LogTrace("PFAlgo|createCandidatesHCAL") << "\t\tsum p = " << totalChargedMomentum << " +- " << sqrt(sumpError2);
2285 LogTrace("PFAlgo|createCandidatesHCAL") << "\t\tsum ecal = " << totalEcal;
2286 LogTrace("PFAlgo|createCandidatesHCAL") << "\t\tsum hcal = " << totalHcal;
2287 LogTrace("PFAlgo|createCandidatesHCAL") << "\t\t => Calo Energy = " << caloEnergy << " +- " << caloResolution;
2288 LogTrace("PFAlgo|createCandidatesHCAL")
2289 << "\t\t => Calo Energy- total charged momentum = " << caloEnergy - totalChargedMomentum << " +- "
2290 << sqrt(sumpError2 + caloResolution * caloResolution);
2291 #endif
2292
2293
2294 unsigned corrTrack = 10000000;
2295 double corrFact = 1.;
2296
2297 if (rejectTracks_Bad_ && totalChargedMomentum - caloEnergy > nSigmaTRACK_ * caloResolution) {
2298 for (auto const& trk : associatedTracks) {
2299 const unsigned iTrack = trk.second.first;
2300
2301 if (!active[iTrack])
2302 continue;
2303 const reco::TrackRef& trackref = elements[trk.second.first].trackRef();
2304
2305 const double dptRel = fabs(trk.first) / trackref->pt() * 100;
2306 const bool isSecondary = isFromSecInt(elements[iTrack], "secondary");
2307 const bool isPrimary = isFromSecInt(elements[iTrack], "primary");
2308
2309 if (isSecondary && dptRel < dptRel_DispVtx_)
2310 continue;
2311
2312 if (fabs(trk.first) < ptError_)
2313 break;
2314
2315 const double wouldBeTotalChargedMomentum = totalChargedMomentum - trackref->p();
2316
2317
2318
2319 if (wouldBeTotalChargedMomentum > caloEnergy) {
2320 if (isSecondary) {
2321 LogTrace("PFAlgo|createCandidatesHCAL")
2322 << "In bad track rejection step dptRel = " << dptRel << " dptRel_DispVtx_ = " << dptRel_DispVtx_;
2323 LogTrace("PFAlgo|createCandidatesHCAL")
2324 << "The calo energy would be still smaller even without this track but it is attached to a NI";
2325 }
2326
2327 if (isPrimary || (isSecondary && dptRel < dptRel_DispVtx_))
2328 continue;
2329 active[iTrack] = false;
2330 totalChargedMomentum = wouldBeTotalChargedMomentum;
2331 LogTrace("PFAlgo|createCandidatesHCAL")
2332 << "\tElement " << elements[iTrack] << " rejected (dpt = " << -trk.first
2333 << " GeV/c, algo = " << trackref->algo() << ")";
2334
2335 } else {
2336 if (isPrimary)
2337 break;
2338 corrTrack = iTrack;
2339 corrFact = (caloEnergy - wouldBeTotalChargedMomentum) / elements[trk.second.first].trackRef()->p();
2340 if (trackref->p() * corrFact < 0.05) {
2341 corrFact = 0.;
2342 active[iTrack] = false;
2343 }
2344 totalChargedMomentum -= trackref->p() * (1. - corrFact);
2345 LogTrace("PFAlgo|createCandidatesHCAL")
2346 << "\tElement " << elements[iTrack] << " (dpt = " << -trk.first << " GeV/c, algo = " << trackref->algo()
2347 << ") rescaled by " << corrFact << " Now the total charged momentum is " << totalChargedMomentum;
2348 break;
2349 }
2350 }
2351 }
2352
2353
2354 caloResolution = neutralHadronEnergyResolution(totalChargedMomentum, hclusterref->positionREP().Eta());
2355 caloResolution *= totalChargedMomentum;
2356 caloResolution = std::sqrt(caloResolution * caloResolution + muonHCALError + muonECALError);
2357
2358
2359
2360
2361 if (rejectTracks_Step45_ && sortedTracks.size() > 1 &&
2362 totalChargedMomentum - caloEnergy > nSigmaTRACK_ * caloResolution) {
2363 for (auto const& trk : associatedTracks) {
2364 unsigned iTrack = trk.second.first;
2365 reco::TrackRef trackref = elements[iTrack].trackRef();
2366 if (!active[iTrack])
2367 continue;
2368
2369 double dptRel = fabs(trk.first) / trackref->pt() * 100;
2370 bool isPrimaryOrSecondary = isFromSecInt(elements[iTrack], "all");
2371
2372 if (isPrimaryOrSecondary && dptRel < dptRel_DispVtx_)
2373 continue;
2374
2375 if (PFTrackAlgoTools::step5(trackref->algo())) {
2376 active[iTrack] = false;
2377 totalChargedMomentum -= trackref->p();
2378
2379 LogTrace("PFAlgo|createCandidatesHCAL")
2380 << "\tElement " << elements[iTrack] << " rejected (dpt = " << -trk.first
2381 << " GeV/c, algo = " << trackref->algo() << ")";
2382 }
2383 }
2384 }
2385
2386
2387 caloResolution = neutralHadronEnergyResolution(totalChargedMomentum, hclusterref->positionREP().Eta());
2388 caloResolution *= totalChargedMomentum;
2389 caloResolution = std::sqrt(caloResolution * caloResolution + muonHCALError + muonECALError);
2390
2391
2392 sumpError2 = 0.;
2393 for (auto const& trk : associatedTracks) {
2394 unsigned iTrack = trk.second.first;
2395 if (!active[iTrack])
2396 continue;
2397 reco::TrackRef trackRef = elements[iTrack].trackRef();
2398 double trackMomentum = trackRef->p();
2399 double dp = trackRef->qoverpError() * trackMomentum * trackMomentum;
2400 unsigned tmpi = reconstructTrack(elements[iTrack]);
2401
2402 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iTrack);
2403 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iHcal);
2404 setHcalDepthInfo((*pfCandidates_)[tmpi], *hclusterref);
2405 auto myEcals = associatedEcals.equal_range(iTrack);
2406 for (auto ii = myEcals.first; ii != myEcals.second; ++ii) {
2407 unsigned iEcal = ii->second.second;
2408 if (active[iEcal])
2409 continue;
2410 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iEcal);
2411 }
2412
2413 if (useHO_) {
2414 auto myHOs = associatedHOs.equal_range(iTrack);
2415 for (auto ii = myHOs.first; ii != myHOs.second; ++ii) {
2416 unsigned iHO = ii->second.second;
2417 if (active[iHO])
2418 continue;
2419 (*pfCandidates_)[tmpi].addElementInBlock(blockref, iHO);
2420 }
2421 }
2422
2423 if (iTrack == corrTrack) {
2424 if (corrFact < 0.)
2425 corrFact = 0.;
2426 auto& candRescale = (*pfCandidates_)[tmpi];
2427 LogTrace("PFAlgo|createCandidatesHCAL")
2428 << "\tBefore rescaling: momentum " << candRescale.p() << " pT " << candRescale.pt() << " energy "
2429 << candRescale.energy() << " mass " << candRescale.mass() << std::endl
2430 << "\tTo rescale by " << corrFact << std::endl;
2431 candRescale.rescaleMomentum(corrFact);
2432 LogTrace("PFAlgo|createCandidatesHCAL")
2433 << "\tRescaled candidate momentum " << candRescale.p() << " pT " << candRescale.pt() << " energy "
2434 << candRescale.energy() << " mass " << candRescale.mass() << std::endl;
2435 trackMomentum *= corrFact;
2436 }
2437 chargedHadronsIndices.push_back(tmpi);
2438 chargedHadronsInBlock.push_back(iTrack);
2439 active[iTrack] = false;
2440 hcalP.push_back(trackMomentum);
2441 hcalDP.push_back(dp);
2442 if (dp / trackMomentum > maxDPovP)
2443 maxDPovP = dp / trackMomentum;
2444 sumpError2 += dp * dp;
2445 }
2446
2447
2448 double totalError = sqrt(sumpError2 + caloResolution * caloResolution);
2449
2450 #ifdef EDM_ML_DEBUG
2451 LogTrace("PFAlgo|createCandidatesHCAL")
2452 << "\tCompare Calo Energy to total charged momentum " << endl
2453 << "\t\tsum p = " << totalChargedMomentum << " +- " << sqrt(sumpError2) << endl
2454 << "\t\tsum ecal = " << totalEcal << endl
2455 << "\t\tsum hcal = " << totalHcal << endl
2456 << "\t\t => Calo Energy = " << caloEnergy << " +- " << caloResolution << endl
2457 << "\t\t => Calo Energy- total charged momentum = " << caloEnergy - totalChargedMomentum << " +- "
2458 << totalError;
2459 #endif
2460
2461
2462
2463
2464 double nsigma = nSigmaHCAL(totalChargedMomentum, hclusterref->positionREP().Eta());
2465
2466 if (abs(totalChargedMomentum - caloEnergy) < nsigma * totalError) {
2467
2468
2469
2470 #ifdef EDM_ML_DEBUG
2471 LogTrace("PFAlgo|createCandidatesHCAL")
2472 << "\t\tcase 1: COMPATIBLE "
2473 << "|Calo Energy- total charged momentum| = " << abs(caloEnergy - totalChargedMomentum) << " < " << nsigma
2474 << " x " << totalError;
2475 if (maxDPovP < 0.1)
2476 LogTrace("PFAlgo|createCandidatesHCAL") << "\t\t\tmax DP/P = " << maxDPovP << " less than 0.1: do nothing ";
2477 else
2478 LogTrace("PFAlgo|createCandidatesHCAL")
2479 << "\t\t\tmax DP/P = " << maxDPovP << " > 0.1: take weighted averages ";
2480 #endif
2481
2482
2483 if (maxDPovP > 0.1) {
2484
2485
2486 int nrows = chargedHadronsIndices.size();
2487 TMatrixTSym<double> a(nrows);
2488 TVectorD b(nrows);
2489 TVectorD check(nrows);
2490 double sigma2E = caloResolution * caloResolution;
2491 for (int i = 0; i < nrows; i++) {
2492 double sigma2i = hcalDP[i] * hcalDP[i];
2493 LogTrace("PFAlgo|createCandidatesHCAL")
2494 << "\t\t\ttrack associated to hcal " << i << " P = " << hcalP[i] << " +- " << hcalDP[i];
2495 a(i, i) = 1. / sigma2i + 1. / sigma2E;
2496 b(i) = hcalP[i] / sigma2i + caloEnergy / sigma2E;
2497 for (int j = 0; j < nrows; j++) {
2498 if (i == j)
2499 continue;
2500 a(i, j) = 1. / sigma2E;
2501 }
2502 }
2503
2504
2505 TDecompChol decomp(a);
2506 bool ok = false;
2507 TVectorD x = decomp.Solve(b, ok);
2508
2509
2510 if (ok) {
2511 for (int i = 0; i < nrows; i++) {
2512
2513 unsigned ich = chargedHadronsIndices[i];
2514 double rescaleFactor = x(i) / hcalP[i];
2515 if (rescaleFactor < 0.)
2516 rescaleFactor = 0.;
2517 auto& candRescale = (*pfCandidates_)[ich];
2518 LogTrace("PFAlgo|createCandidatesHCAL")
2519 << "\tBefore rescaling: momentum " << candRescale.p() << " pT " << candRescale.pt() << " energy "
2520 << candRescale.energy() << " mass " << candRescale.mass() << std::endl
2521 << "\tTo rescale by " << rescaleFactor << std::endl;
2522 candRescale.rescaleMomentum(rescaleFactor);
2523 LogTrace("PFAlgo|createCandidatesHCAL")
2524 << "\tRescaled candidate momentum " << candRescale.p() << " pT " << candRescale.pt() << " energy "
2525 << candRescale.energy() << " mass " << candRescale.mass() << std::endl;
2526
2527 LogTrace("PFAlgo|createCandidatesHCAL")
2528 << "\t\t\told p " << hcalP[i] << " new p " << x(i) << " rescale " << rescaleFactor;
2529 }
2530 } else {
2531 edm::LogError("PFAlgo::createCandidatesHCAL") << "TDecompChol.Solve returned ok=false";
2532 assert(0);
2533 }
2534 }
2535 }
2536
2537
2538 else if (caloEnergy > totalChargedMomentum) {
2539
2540
2541
2542
2543 double eNeutralHadron = caloEnergy - totalChargedMomentum;
2544 double ePhoton = (caloEnergy - totalChargedMomentum) /
2545 slopeEcal;
2546
2547
2548
2549 #ifdef EDM_ML_DEBUG
2550 if (!sortedTracks.empty()) {
2551 LogTrace("PFAlgo|createCandidatesHCAL") << "\t\tcase 2: NEUTRAL DETECTION " << caloEnergy << " > " << nsigma
2552 << "x" << totalError << " + " << totalChargedMomentum;
2553 LogTrace("PFAlgo|createCandidatesHCAL") << "\t\tneutral activity detected: " << endl
2554 << "\t\t\t photon = " << ePhoton << endl
2555 << "\t\t\tor neutral hadron = " << eNeutralHadron;
2556
2557 LogTrace("PFAlgo|createCandidatesHCAL") << "\t\tphoton or hadron ?";
2558 }
2559
2560 if (sortedTracks.empty())
2561 LogTrace("PFAlgo|createCandidatesHCAL") << "\t\tno track -> hadron ";
2562 else
2563 LogTrace("PFAlgo|createCandidatesHCAL")
2564 << "\t\t" << sortedTracks.size() << " tracks -> check if the excess is photonic or hadronic";
2565 #endif
2566
2567 double ratioMax = 0.;
2568 reco::PFClusterRef maxEcalRef;
2569 unsigned maxiEcal = 9999;
2570
2571
2572
2573 LogTrace("PFAlgo|createCandidatesHCAL") << "loop over sortedTracks.size()=" << sortedTracks.size();
2574 for (auto const& trk : sortedTracks) {
2575 unsigned iTrack = trk.second;
2576
2577 PFBlockElement::Type type = elements[iTrack].type();
2578 assert(type == reco::PFBlockElement::TRACK);
2579
2580 reco::TrackRef trackRef = elements[iTrack].trackRef();
2581 assert(!trackRef.isNull());
2582
2583 auto iae = associatedEcals.find(iTrack);
2584
2585 if (iae == associatedEcals.end())
2586 continue;
2587
2588
2589 unsigned iEcal = iae->second.second;
2590
2591 PFBlockElement::Type typeEcal = elements[iEcal].type();
2592 assert(typeEcal == reco::PFBlockElement::ECAL);
2593
2594 reco::PFClusterRef clusterRef = elements[iEcal].clusterRef();
2595 assert(!clusterRef.isNull());
2596
2597 double pTrack = trackRef->p();
2598 double eECAL = clusterRef->energy();
2599 double eECALOverpTrack = eECAL / pTrack;
2600
2601 if (eECALOverpTrack > ratioMax) {
2602 ratioMax = eECALOverpTrack;
2603 maxEcalRef = clusterRef;
2604 maxiEcal = iEcal;
2605 }
2606
2607 }
2608
2609 std::vector<reco::PFClusterRef> pivotalClusterRef;
2610 std::vector<unsigned> iPivotal;
2611 std::vector<double> particleEnergy, ecalEnergy, hcalEnergy, rawecalEnergy, rawhcalEnergy;
2612 std::vector<::math::XYZVector> particleDirection;
2613
2614
2615
2616 if (ePhoton < totalEcal || eNeutralHadron - calibEcal < 1E-10) {
2617 if (!maxEcalRef.isNull()) {
2618
2619 mergedPhotonEnergy = ePhoton;
2620 }
2621 } else {
2622
2623 if (!maxEcalRef.isNull()) {
2624
2625 mergedPhotonEnergy = totalEcalEGMCalib;
2626 }
2627
2628 mergedNeutralHadronEnergy = eNeutralHadron - calibEcal;
2629 }
2630
2631 if (mergedPhotonEnergy > 0) {
2632
2633
2634
2635 if (ecalClusters.size() <= 1) {
2636 ecalClusters.clear();
2637 ecalClusters.emplace_back(
2638 maxiEcal,
2639 photonAtECAL,
2640 1.);
2641 sumEcalClusters = sqrt(photonAtECAL.Mag2());
2642 }
2643 for (auto const& pae : ecalClusters) {
2644 const double clusterEnergyCalibrated =
2645 sqrt(std::get<1>(pae).Mag2()) *
2646 std::get<2>(
2647 pae);
2648 particleEnergy.push_back(mergedPhotonEnergy * clusterEnergyCalibrated / sumEcalClusters);
2649 particleDirection.push_back(std::get<1>(pae));
2650 ecalEnergy.push_back(mergedPhotonEnergy * clusterEnergyCalibrated / sumEcalClusters);
2651 hcalEnergy.push_back(0.);
2652 rawecalEnergy.push_back(totalEcal);
2653 rawhcalEnergy.push_back(0.);
2654 pivotalClusterRef.push_back(elements[std::get<0>(pae)].clusterRef());
2655 iPivotal.push_back(std::get<0>(pae));
2656 }
2657 }
2658
2659 if (mergedNeutralHadronEnergy > 1.0) {
2660
2661
2662 if (ecalClusters.size() <= 1) {
2663 ecalClusters.clear();
2664 ecalClusters.emplace_back(
2665 iHcal,
2666 hadronAtECAL,
2667 1.);
2668 sumEcalClusters = sqrt(hadronAtECAL.Mag2());
2669 }
2670 for (auto const& pae : ecalClusters) {
2671 const double clusterEnergyCalibrated =
2672 sqrt(std::get<1>(pae).Mag2()) *
2673 std::get<2>(
2674 pae);
2675 particleEnergy.push_back(mergedNeutralHadronEnergy * clusterEnergyCalibrated / sumEcalClusters);
2676 particleDirection.push_back(std::get<1>(pae));
2677 ecalEnergy.push_back(0.);
2678 hcalEnergy.push_back(mergedNeutralHadronEnergy * clusterEnergyCalibrated / sumEcalClusters);
2679 rawecalEnergy.push_back(0.);
2680 rawhcalEnergy.push_back(totalHcal);
2681 pivotalClusterRef.push_back(hclusterref);
2682 iPivotal.push_back(iHcal);
2683 }
2684 }
2685
2686
2687
2688
2689
2690 for (unsigned iPivot = 0; iPivot < iPivotal.size(); ++iPivot) {
2691 if (particleEnergy[iPivot] < 0.)
2692 edm::LogWarning("PFAlgo|createCandidatesHCAL")
2693 << "ALARM = Negative energy for iPivot=" << iPivot << ", " << particleEnergy[iPivot];
2694
2695 const bool useDirection = true;
2696 auto& neutral = (*pfCandidates_)[reconstructCluster(*pivotalClusterRef[iPivot],
2697 particleEnergy[iPivot],
2698 useDirection,
2699 particleDirection[iPivot].X(),
2700 particleDirection[iPivot].Y(),
2701 particleDirection[iPivot].Z())];
2702
2703 neutral.setEcalEnergy(rawecalEnergy[iPivot], ecalEnergy[iPivot]);
2704 if (!useHO_) {
2705 neutral.setHcalEnergy(rawhcalEnergy[iPivot], hcalEnergy[iPivot]);
2706 neutral.setHoEnergy(0., 0.);
2707 } else {
2708 if (rawhcalEnergy[iPivot] == 0.) {
2709 neutral.setHcalEnergy(0., 0.);
2710 neutral.setHoEnergy(0., 0.);
2711 } else {
2712 neutral.setHcalEnergy(max(rawhcalEnergy[iPivot] - totalHO, 0.0),
2713 hcalEnergy[iPivot] * max(1. - totalHO / rawhcalEnergy[iPivot], 0.));
2714 neutral.setHoEnergy(totalHO, totalHO * hcalEnergy[iPivot] / rawhcalEnergy[iPivot]);
2715 }
2716 }
2717 neutral.setPs1Energy(0.);
2718 neutral.setPs2Energy(0.);
2719 neutral.set_mva_nothing_gamma(-1.);
2720
2721
2722 neutral.addElementInBlock(blockref, iHcal);
2723 for (unsigned iTrack : chargedHadronsInBlock) {
2724 neutral.addElementInBlock(blockref, iTrack);
2725
2726 const ::math::XYZPointF& chargedPosition =
2727 dynamic_cast<const reco::PFBlockElementTrack*>(&elements[iTrack])->positionAtECALEntrance();
2728 neutral.setPositionAtECALEntrance(chargedPosition);
2729
2730 auto myEcals = associatedEcals.equal_range(iTrack);
2731 for (auto ii = myEcals.first; ii != myEcals.second; ++ii) {
2732 unsigned iEcal = ii->second.second;
2733 if (active[iEcal])
2734 continue;
2735 neutral.addElementInBlock(blockref, iEcal);
2736 }
2737 }
2738 }
2739
2740 }
2741
2742
2743
2744
2745
2746
2747
2748
2749
2750 double totalHcalEnergyCalibrated = std::max(calibHcal - mergedNeutralHadronEnergy, 0.);
2751
2752 double totalEcalEnergyCalibrated = std::max(calibEcal - mergedPhotonEnergy, 0.);
2753
2754
2755
2756
2757 double chargedHadronsTotalEnergy = 0;
2758 for (unsigned index : chargedHadronsIndices) {
2759 reco::PFCandidate& chargedHadron = (*pfCandidates_)[index];
2760 chargedHadronsTotalEnergy += chargedHadron.energy();
2761 }
2762
2763 for (unsigned index : chargedHadronsIndices) {
2764 reco::PFCandidate& chargedHadron = (*pfCandidates_)[index];
2765 float fraction = chargedHadron.energy() / chargedHadronsTotalEnergy;
2766
2767 if (!useHO_) {
2768 chargedHadron.setHcalEnergy(fraction * totalHcal, fraction * totalHcalEnergyCalibrated);
2769 chargedHadron.setHoEnergy(0., 0.);
2770 } else {
2771 chargedHadron.setHcalEnergy(fraction * max(totalHcal - totalHO, 0.0),
2772 fraction * totalHcalEnergyCalibrated * (1. - totalHO / totalHcal));
2773 chargedHadron.setHoEnergy(fraction * totalHO, fraction * totalHO * totalHcalEnergyCalibrated / totalHcal);
2774 }
2775
2776 chargedHadron.setEcalEnergy(fraction * totalEcal, fraction * totalEcalEnergyCalibrated);
2777 }
2778
2779
2780 for (auto const& ecalSatellite : ecalSatellites) {
2781
2782 unsigned iEcal = std::get<0>(ecalSatellite.second);
2783 if (!active[iEcal])
2784 continue;
2785
2786
2787 PFBlockElement::Type type = elements[iEcal].type();
2788 assert(type == PFBlockElement::ECAL);
2789 PFClusterRef eclusterref = elements[iEcal].clusterRef();
2790 assert(!eclusterref.isNull());
2791
2792
2793 active[iEcal] = false;
2794
2795
2796 std::multimap<double, unsigned> assTracks;
2797 block.associatedElements(iEcal, linkData, assTracks, reco::PFBlockElement::TRACK, reco::PFBlock::LINKTEST_ALL);
2798
2799
2800 double ecalClusterEnergyCalibrated =
2801 sqrt(std::get<1>(ecalSatellite.second).Mag2()) *
2802 std::get<2>(
2803 ecalSatellite.second);
2804 auto& cand = (*pfCandidates_)[reconstructCluster(*eclusterref, ecalClusterEnergyCalibrated)];
2805 cand.setEcalEnergy(eclusterref->energy(), ecalClusterEnergyCalibrated);
2806 cand.setHcalEnergy(0., 0.);
2807 cand.setHoEnergy(0., 0.);
2808 cand.setPs1Energy(associatedPSs[iEcal].first);
2809 cand.setPs2Energy(associatedPSs[iEcal].second);
2810 cand.addElementInBlock(blockref, iEcal);
2811 cand.addElementInBlock(blockref, sortedTracks.begin()->second);
2812
2813 if (fabs(eclusterref->energy() - sqrt(std::get<1>(ecalSatellite.second).Mag2())) > 1e-3 ||
2814 fabs(eclusterref->correctedEnergy() - ecalClusterEnergyCalibrated) > 1e-3)
2815 edm::LogWarning("PFAlgo:processBlock")
2816 << "ecalCluster vs ecalSatellites look inconsistent (eCluster E, calibE, ecalSatellite E, calib E): "
2817 << eclusterref->energy() << " " << eclusterref->correctedEnergy() << " "
2818 << sqrt(std::get<1>(ecalSatellite.second).Mag2()) << " " << ecalClusterEnergyCalibrated;
2819
2820 }
2821
2822 }
2823
2824 LogTrace("PFAlgo|createCandidatesHCAL") << "end of function PFAlgo::createCandidatesHCAL";
2825 }
2826
2827 void PFAlgo::createCandidatesHCALUnlinked(const reco::PFBlock& block,
2828 reco::PFBlock::LinkData& linkData,
2829 const edm::OwnVector<reco::PFBlockElement>& elements,
2830 std::vector<bool>& active,
2831 const reco::PFBlockRef& blockref,
2832 ElementIndices& inds,
2833 std::vector<bool>& deadArea) {
2834
2835 LogTrace("PFAlgo|createCandidatesHCALUnlinked")
2836 << "start of function PFAlgo::createCandidatesHCALUnlinked, hcalIs.size()=" << inds.hcalIs.size();
2837
2838
2839
2840 for (unsigned iHcal : inds.hcalIs) {
2841
2842 std::vector<unsigned> ecalRefs;
2843 std::vector<unsigned> hoRefs;
2844
2845 LogTrace("PFAlgo|createCandidatesHCALUnlinked") << elements[iHcal] << " ";
2846
2847 if (!active[iHcal]) {
2848 LogTrace("PFAlgo|createCandidatesHCALUnlinked") << "not active " << iHcal;
2849 continue;
2850 }
2851
2852
2853 std::multimap<double, unsigned> ecalElems;
2854 block.associatedElements(iHcal, linkData, ecalElems, reco::PFBlockElement::ECAL, reco::PFBlock::LINKTEST_ALL);
2855
2856
2857 float totalEcal = 0.;
2858 float ecalMax = 0.;
2859 reco::PFClusterRef eClusterRef;
2860 for (auto const& ecal : ecalElems) {
2861 unsigned iEcal = ecal.second;
2862 double dist = ecal.first;
2863 PFBlockElement::Type type = elements[iEcal].type();
2864 assert(type == PFBlockElement::ECAL);
2865
2866
2867 if (!active[iEcal])
2868 continue;
2869
2870
2871
2872
2873
2874
2875
2876
2877
2878
2879
2880
2881
2882
2883
2884
2885
2886
2887
2888
2889 std::multimap<double, unsigned> hcalElems;
2890 block.associatedElements(iEcal, linkData, hcalElems, reco::PFBlockElement::HCAL, reco::PFBlock::LINKTEST_ALL);
2891
2892 const bool isClosest = std::none_of(hcalElems.begin(), hcalElems.end(), [&](auto const& hcal) {
2893 return active[hcal.second] && hcal.first < dist;
2894 });
2895
2896 if (!isClosest)
2897 continue;
2898
2899 #ifdef EDM_ML_DEBUG
2900 LogTrace("PFAlgo|createCandidatesHCALUnlinked")
2901 << "\telement " << elements[iEcal] << " linked with dist " << dist;
2902 LogTrace("PFAlgo|createCandidatesHCALUnlinked") << "Added to HCAL cluster to form a neutral hadron";
2903 #endif
2904
2905 reco::PFClusterRef eclusterRef = elements[iEcal].clusterRef();
2906 assert(!eclusterRef.isNull());
2907
2908
2909 double ecalEnergy = eclusterRef->energy();
2910
2911 totalEcal += ecalEnergy;
2912 if (ecalEnergy > ecalMax) {
2913 ecalMax = ecalEnergy;
2914 eClusterRef = eclusterRef;
2915 }
2916
2917 ecalRefs.push_back(iEcal);
2918 active[iEcal] = false;
2919
2920 }
2921
2922
2923 double totalHO = 0.;
2924 double hoMax = 0.;
2925
2926 if (useHO_) {
2927 std::multimap<double, unsigned> hoElems;
2928 block.associatedElements(iHcal, linkData, hoElems, reco::PFBlockElement::HO, reco::PFBlock::LINKTEST_ALL);
2929
2930
2931
2932
2933
2934 reco::PFClusterRef hoClusterRef;
2935 for (auto const& ho : hoElems) {
2936 unsigned iHO = ho.second;
2937 double dist = ho.first;
2938 PFBlockElement::Type type = elements[iHO].type();
2939 assert(type == PFBlockElement::HO);
2940
2941
2942 if (!active[iHO])
2943 continue;
2944
2945
2946
2947
2948
2949 std::multimap<double, unsigned> hcalElems;
2950 block.associatedElements(iHO, linkData, hcalElems, reco::PFBlockElement::HCAL, reco::PFBlock::LINKTEST_ALL);
2951
2952 const bool isClosest = std::none_of(hcalElems.begin(), hcalElems.end(), [&](auto const& hcal) {
2953 return active[hcal.second] && hcal.first < dist;
2954 });
2955
2956 if (!isClosest)
2957 continue;
2958
2959 #ifdef EDM_ML_DEBUG
2960 if (useHO_) {
2961 LogTrace("PFAlgo|createCandidatesHCALUnlinked")
2962 << "\telement " << elements[iHO] << " linked with dist " << dist;
2963 LogTrace("PFAlgo|createCandidatesHCALUnlinked") << "Added to HCAL cluster to form a neutral hadron";
2964 }
2965 #endif
2966
2967 reco::PFClusterRef hoclusterRef = elements[iHO].clusterRef();
2968 assert(!hoclusterRef.isNull());
2969
2970 double hoEnergy =
2971 hoclusterRef->energy();
2972
2973 totalHO += hoEnergy;
2974 if (hoEnergy > hoMax) {
2975 hoMax = hoEnergy;
2976 hoClusterRef = hoclusterRef;
2977
2978 }
2979
2980 hoRefs.push_back(iHO);
2981 active[iHO] = false;
2982
2983 }
2984 }
2985
2986 PFClusterRef hclusterRef = elements[iHcal].clusterRef();
2987 assert(!hclusterRef.isNull());
2988
2989
2990 double totalHcal = hclusterRef->energy();
2991
2992 if (useHO_)
2993 totalHcal += totalHO;
2994
2995
2996 double calibEcal = totalEcal > 0. ? totalEcal : 0.;
2997 double calibHcal = std::max(0., totalHcal);
2998 if (hclusterRef->layer() == PFLayer::HF_HAD || hclusterRef->layer() == PFLayer::HF_EM) {
2999 calibEcal = totalEcal;
3000 } else {
3001 calibration_.energyEmHad(
3002 -1., calibEcal, calibHcal, hclusterRef->positionREP().Eta(), hclusterRef->positionREP().Phi());
3003 }
3004
3005 auto& cand = (*pfCandidates_)[reconstructCluster(*hclusterRef, calibEcal + calibHcal)];
3006
3007 cand.setEcalEnergy(totalEcal, calibEcal);
3008 if (!useHO_) {
3009 cand.setHcalEnergy(totalHcal, calibHcal);
3010 cand.setHoEnergy(0., 0.);
3011 } else {
3012 cand.setHcalEnergy(max(totalHcal - totalHO, 0.0), calibHcal * (1. - totalHO / totalHcal));
3013 cand.setHoEnergy(totalHO, totalHO * calibHcal / totalHcal);
3014 }
3015 cand.setPs1Energy(0.);
3016 cand.setPs2Energy(0.);
3017 cand.addElementInBlock(blockref, iHcal);
3018 for (auto const& ec : ecalRefs)
3019 cand.addElementInBlock(blockref, ec);
3020 for (auto const& ho : hoRefs)
3021 cand.addElementInBlock(blockref, ho);
3022
3023 }
3024 }
3025
3026 void PFAlgo::createCandidatesECAL(const reco::PFBlock& block,
3027 reco::PFBlock::LinkData& linkData,
3028 const edm::OwnVector<reco::PFBlockElement>& elements,
3029 std::vector<bool>& active,
3030 const reco::PFBlockRef& blockref,
3031 ElementIndices& inds,
3032 std::vector<bool>& deadArea) {
3033 LogTrace("PFAlgo|createCandidatesECAL")
3034 << "start of function PFAlgo::createCandidatesECAL(), ecalIs.size()=" << inds.ecalIs.size();
3035
3036
3037
3038
3039
3040 for (unsigned i = 0; i < inds.ecalIs.size(); i++) {
3041 unsigned iEcal = inds.ecalIs[i];
3042
3043 LogTrace("PFAlgo|createCandidatesECAL") << "elements[" << iEcal << "]=" << elements[iEcal];
3044
3045 if (!active[iEcal]) {
3046 LogTrace("PFAlgo|createCandidatesECAL") << "iEcal=" << iEcal << " not active";
3047 continue;
3048 }
3049
3050 PFBlockElement::Type type = elements[iEcal].type();
3051 assert(type == PFBlockElement::ECAL);
3052
3053 PFClusterRef clusterref = elements[iEcal].clusterRef();
3054 assert(!clusterref.isNull());
3055
3056 active[iEcal] = false;
3057
3058 float ecalEnergy = clusterref->correctedEnergy();
3059
3060 double particleEnergy = ecalEnergy;
3061
3062 auto& cand = (*pfCandidates_)[reconstructCluster(*clusterref, particleEnergy)];
3063
3064 cand.setEcalEnergy(clusterref->energy(), ecalEnergy);
3065 cand.setHcalEnergy(0., 0.);
3066 cand.setHoEnergy(0., 0.);
3067 cand.setPs1Energy(0.);
3068 cand.setPs2Energy(0.);
3069 cand.addElementInBlock(blockref, iEcal);
3070
3071 }
3072 LogTrace("PFAlgo|createCandidatesECAL") << "end of function PFALgo::createCandidatesECAL";
3073 }
3074
3075 void PFAlgo::processBlock(const reco::PFBlockRef& blockref,
3076 std::list<reco::PFBlockRef>& hcalBlockRefs,
3077 std::list<reco::PFBlockRef>& ecalBlockRefs,
3078 PFEGammaFilters const* pfegamma) {
3079 assert(!blockref.isNull());
3080 const reco::PFBlock& block = *blockref;
3081
3082 LogTrace("PFAlgo|processBlock") << "start of function PFAlgo::processBlock, block=" << block;
3083
3084 const edm::OwnVector<reco::PFBlockElement>& elements = block.elements();
3085 LogTrace("PFAlgo|processBlock") << "elements.size()=" << elements.size();
3086
3087 PFBlock::LinkData linkData = block.linkData();
3088
3089
3090 vector<bool> active(elements.size(), true);
3091
3092
3093
3094 std::vector<reco::PFCandidate> tempElectronCandidates;
3095 tempElectronCandidates.clear();
3096
3097
3098 if (useEGammaFilters_) {
3099 egammaFilters(blockref, active, pfegamma);
3100 }
3101
3102
3103 if (usePFConversions_) {
3104 conversionAlgo(elements, active);
3105 }
3106
3107
3108
3109
3110
3111
3112
3113
3114
3115
3116
3117
3118
3119
3120
3121
3122
3123
3124
3125
3126
3127
3128 vector<bool> deadArea(elements.size(), false);
3129
3130
3131 ElementIndices inds;
3132
3133 elementLoop(block, linkData, elements, active, blockref, inds, deadArea);
3134
3135
3136
3137 if (!(inds.hfEmIs.empty() && inds.hfHadIs.empty())) {
3138 createCandidatesHF(block, linkData, elements, active, blockref, inds);
3139 if (inds.hcalIs.empty() && inds.ecalIs.empty())
3140 return;
3141 LogDebug("PFAlgo::processBlock")
3142 << "Block contains HF clusters, and also contains ECAL or HCAL clusters. Continue.\n"
3143 << block;
3144 }
3145
3146 createCandidatesHCAL(block, linkData, elements, active, blockref, inds, deadArea);
3147
3148 createCandidatesHCALUnlinked(block, linkData, elements, active, blockref, inds, deadArea);
3149 createCandidatesECAL(block, linkData, elements, active, blockref, inds, deadArea);
3150
3151 LogTrace("PFAlgo|processBlock") << "end of function PFAlgo::processBlock";
3152 }
3153
3154
3155 unsigned PFAlgo::reconstructTrack(const reco::PFBlockElement& elt, bool allowLoose) {
3156 const auto* eltTrack = dynamic_cast<const reco::PFBlockElementTrack*>(&elt);
3157
3158 const reco::TrackRef& trackRef = eltTrack->trackRef();
3159 const reco::Track& track = *trackRef;
3160 const reco::MuonRef& muonRef = eltTrack->muonRef();
3161 int charge = track.charge() > 0 ? 1 : -1;
3162
3163
3164 double px = track.px();
3165 double py = track.py();
3166 double pz = track.pz();
3167 double energy = sqrt(track.p() * track.p() + 0.13957 * 0.13957);
3168
3169 LogTrace("PFAlgo|reconstructTrack") << "Reconstructing PF candidate from track of pT = " << track.pt()
3170 << " eta = " << track.eta() << " phi = " << track.phi() << " px = " << px
3171 << " py = " << py << " pz = " << pz << " energy = " << energy;
3172
3173
3174 ::math::XYZTLorentzVector momentum(px, py, pz, energy);
3175 reco::PFCandidate::ParticleType particleType = reco::PFCandidate::h;
3176
3177
3178 LogTrace("PFAlgo|reconstructTrack") << "Creating PFCandidate charge=" << charge << ", type=" << particleType
3179 << ", pt=" << momentum.pt() << ", eta=" << momentum.eta()
3180 << ", phi=" << momentum.phi();
3181 pfCandidates_->push_back(PFCandidate(charge, momentum, particleType));
3182
3183 pfCandidates_->back().setVertex(trackRef->vertex());
3184 pfCandidates_->back().setTrackRef(trackRef);
3185 pfCandidates_->back().setPositionAtECALEntrance(eltTrack->positionAtECALEntrance());
3186 if (muonRef.isNonnull())
3187 pfCandidates_->back().setMuonRef(muonRef);
3188
3189
3190 if (elt.isTimeValid())
3191 pfCandidates_->back().setTime(elt.time(), elt.timeError());
3192
3193
3194 bool isMuon = pfmu_->reconstructMuon(pfCandidates_->back(), muonRef, allowLoose);
3195 bool isFromDisp = isFromSecInt(elt, "secondary");
3196
3197 if ((!isMuon) && isFromDisp) {
3198 double dpt = trackRef->ptError();
3199 double dptRel = dpt / trackRef->pt() * 100;
3200
3201
3202
3203 if (dptRel < dptRel_DispVtx_) {
3204 LogTrace("PFAlgo|reconstructTrack")
3205 << "Not refitted px = " << px << " py = " << py << " pz = " << pz << " energy = " << energy;
3206
3207 reco::PFDisplacedVertexRef vRef =
3208 eltTrack->displacedVertexRef(reco::PFBlockElement::T_FROM_DISP)->displacedVertexRef();
3209 reco::Track trackRefit = vRef->refittedTrack(trackRef);
3210
3211 ::math::XYZTLorentzVector momentum(
3212 trackRefit.px(), trackRefit.py(), trackRefit.pz(), sqrt(trackRefit.p() * trackRefit.p() + 0.13957 * 0.13957));
3213 LogTrace("PFAlgo|reconstructTrack")
3214 << "Refitted px = " << px << " py = " << py << " pz = " << pz << " energy = " << energy;
3215 }
3216 pfCandidates_->back().setFlag(reco::PFCandidate::T_FROM_DISP, true);
3217 pfCandidates_->back().setDisplacedVertexRef(
3218 eltTrack->displacedVertexRef(reco::PFBlockElement::T_FROM_DISP)->displacedVertexRef(),
3219 reco::PFCandidate::T_FROM_DISP);
3220 }
3221
3222
3223 if (isFromSecInt(elt, "primary") && !isMuon) {
3224 pfCandidates_->back().setFlag(reco::PFCandidate::T_TO_DISP, true);
3225 pfCandidates_->back().setDisplacedVertexRef(
3226 eltTrack->displacedVertexRef(reco::PFBlockElement::T_TO_DISP)->displacedVertexRef(),
3227 reco::PFCandidate::T_TO_DISP);
3228 }
3229
3230
3231 return pfCandidates_->size() - 1;
3232 }
3233
3234 unsigned PFAlgo::reconstructCluster(const reco::PFCluster& cluster,
3235 double particleEnergy,
3236 bool useDirection,
3237 double particleX,
3238 double particleY,
3239 double particleZ) {
3240 LogTrace("PFAlgo|reconstructCluster") << "start of function PFAlgo::reconstructCluster, cluster=" << cluster
3241 << "particleEnergy=" << particleEnergy << "useDirection=" << useDirection
3242 << "particleX=" << particleX << "particleY=" << particleY
3243 << "particleZ=" << particleZ;
3244
3245 reco::PFCandidate::ParticleType particleType = reco::PFCandidate::X;
3246
3247
3248
3249
3250
3251 double factor = 1.;
3252 if (useDirection) {
3253 switch (cluster.layer()) {
3254 case PFLayer::ECAL_BARREL:
3255 case PFLayer::HCAL_BARREL1:
3256 factor = std::sqrt(cluster.position().Perp2() / (particleX * particleX + particleY * particleY));
3257 break;
3258 case PFLayer::ECAL_ENDCAP:
3259 case PFLayer::HCAL_ENDCAP:
3260 case PFLayer::HF_HAD:
3261 case PFLayer::HF_EM:
3262 factor = cluster.position().Z() / particleZ;
3263 break;
3264 default:
3265 assert(0);
3266 }
3267 }
3268
3269 ::math::XYZPoint vertexPos;
3270 if (useVertices_)
3271 vertexPos = ::math::XYZPoint(primaryVertex_.x(), primaryVertex_.y(), primaryVertex_.z());
3272 else
3273 vertexPos = ::math::XYZPoint(0.0, 0.0, 0.0);
3274
3275 ::math::XYZVector clusterPos(cluster.position().X() - vertexPos.X(),
3276 cluster.position().Y() - vertexPos.Y(),
3277 cluster.position().Z() - vertexPos.Z());
3278 ::math::XYZVector particleDirection(
3279 particleX * factor - vertexPos.X(), particleY * factor - vertexPos.Y(), particleZ * factor - vertexPos.Z());
3280
3281
3282
3283
3284 clusterPos = useDirection ? particleDirection.Unit() : clusterPos.Unit();
3285 clusterPos *= particleEnergy;
3286
3287
3288
3289
3290 double mass = 0;
3291 ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double>> momentum(
3292 clusterPos.X(), clusterPos.Y(), clusterPos.Z(), mass);
3293
3294 ::math::XYZTLorentzVector tmp;
3295
3296 tmp = momentum;
3297
3298
3299 int charge = 0;
3300
3301
3302 switch (cluster.layer()) {
3303 case PFLayer::ECAL_BARREL:
3304 case PFLayer::ECAL_ENDCAP:
3305 particleType = PFCandidate::gamma;
3306 break;
3307 case PFLayer::HCAL_BARREL1:
3308 case PFLayer::HCAL_ENDCAP:
3309 particleType = PFCandidate::h0;
3310 break;
3311 case PFLayer::HF_HAD:
3312 particleType = PFCandidate::h_HF;
3313 break;
3314 case PFLayer::HF_EM:
3315 particleType = PFCandidate::egamma_HF;
3316 break;
3317 default:
3318 assert(0);
3319 }
3320
3321
3322 LogTrace("PFAlgo|reconstructCluster") << "Creating PFCandidate charge=" << charge << ", type=" << particleType
3323 << ", pt=" << tmp.pt() << ", eta=" << tmp.eta() << ", phi=" << tmp.phi();
3324 pfCandidates_->push_back(PFCandidate(charge, tmp, particleType));
3325
3326
3327
3328 pfCandidates_->back().setPositionAtECALEntrance(
3329 ::math::XYZPointF(cluster.position().X(), cluster.position().Y(), cluster.position().Z()));
3330
3331
3332 pfCandidates_->back().setVertex(vertexPos);
3333
3334
3335 setHcalDepthInfo(pfCandidates_->back(), cluster);
3336
3337
3338
3339 LogTrace("PFAlgo|reconstructCluster") << "** candidate: " << pfCandidates_->back();
3340
3341
3342 return pfCandidates_->size() - 1;
3343 }
3344
3345 void PFAlgo::setHcalDepthInfo(reco::PFCandidate& cand, const reco::PFCluster& cluster) const {
3346 std::array<double, 7> energyPerDepth;
3347 std::fill(energyPerDepth.begin(), energyPerDepth.end(), 0.0);
3348 for (auto& hitRefAndFrac : cluster.recHitFractions()) {
3349 const auto& hit = *hitRefAndFrac.recHitRef();
3350 if (DetId(hit.detId()).det() == DetId::Hcal) {
3351 if (hit.depth() == 0) {
3352 edm::LogWarning("setHcalDepthInfo") << "Depth zero found";
3353 continue;
3354 }
3355 if (hit.depth() < 1 || hit.depth() > 7) {
3356 throw cms::Exception("CorruptData") << "Bogus depth " << hit.depth() << " at detid " << hit.detId() << "\n";
3357 }
3358 energyPerDepth[hit.depth() - 1] += hitRefAndFrac.fraction() * hit.energy();
3359 }
3360 }
3361 double sum = std::accumulate(energyPerDepth.begin(), energyPerDepth.end(), 0.);
3362 std::array<float, 7> depthFractions;
3363 if (sum > 0) {
3364 for (unsigned int i = 0; i < depthFractions.size(); ++i) {
3365 depthFractions[i] = energyPerDepth[i] / sum;
3366 }
3367 } else {
3368 std::fill(depthFractions.begin(), depthFractions.end(), 0.f);
3369 }
3370 cand.setHcalDepthEnergyFractions(depthFractions);
3371 }
3372
3373
3374
3375 double PFAlgo::neutralHadronEnergyResolution(double clusterEnergyHCAL, double eta) const {
3376
3377 clusterEnergyHCAL = std::max(clusterEnergyHCAL, 1.);
3378
3379 double resol = fabs(eta) < 1.48 ? sqrt(1.02 * 1.02 / clusterEnergyHCAL + 0.065 * 0.065)
3380 : sqrt(1.20 * 1.20 / clusterEnergyHCAL + 0.028 * 0.028);
3381
3382 return resol;
3383 }
3384
3385 double PFAlgo::nSigmaHCAL(double clusterEnergyHCAL, double eta) const {
3386 double nS = fabs(eta) < 1.48 ? nSigmaHCAL_ * (1. + exp(-clusterEnergyHCAL / nSigmaEConstHCAL))
3387 : nSigmaHCAL_ * (1. + exp(-clusterEnergyHCAL / nSigmaEConstHCAL));
3388
3389 return nS;
3390 }
3391
3392 double PFAlgo::hfEnergyResolution(double clusterEnergyHF) const {
3393
3394 clusterEnergyHF = std::max(clusterEnergyHF, 1.);
3395
3396 double resol =
3397 sqrt(resolHF_square_[0] / clusterEnergyHF + resolHF_square_[1] + resolHF_square_[2] / pow(clusterEnergyHF, 2));
3398
3399
3400
3401 return resol;
3402 }
3403
3404 double PFAlgo::nSigmaHFEM(double clusterEnergyHF) const {
3405 double nS = nSigmaHFEM_ * (1. + exp(-clusterEnergyHF / nSigmaEConstHFEM));
3406 return nS;
3407 }
3408
3409 double PFAlgo::nSigmaHFHAD(double clusterEnergyHF) const {
3410 double nS = nSigmaHFHAD_ * (1. + exp(-clusterEnergyHF / nSigmaEConstHFHAD));
3411 return nS;
3412 }
3413
3414 ostream& operator<<(ostream& out, const PFAlgo& algo) {
3415 if (!out)
3416 return out;
3417
3418 out << "====== Particle Flow Algorithm ======= ";
3419 out << endl;
3420 out << "nSigmaECAL_ " << algo.nSigmaECAL_ << endl;
3421 out << "nSigmaHCAL_ " << algo.nSigmaHCAL_ << endl;
3422 out << "nSigmaHFEM_ " << algo.nSigmaHFEM_ << endl;
3423 out << "nSigmaHFHAD_ " << algo.nSigmaHFHAD_ << endl;
3424 out << endl;
3425 out << algo.calibration_ << endl;
3426 out << endl;
3427 out << "reconstructed particles: " << endl;
3428
3429 if (!algo.pfCandidates_.get()) {
3430 out << "candidates already transfered" << endl;
3431 return out;
3432 }
3433 for (auto const& c : *algo.pfCandidates_)
3434 out << c << endl;
3435
3436 return out;
3437 }
3438
3439 void PFAlgo::associatePSClusters(unsigned iEcal,
3440 reco::PFBlockElement::Type psElementType,
3441 const reco::PFBlock& block,
3442 const edm::OwnVector<reco::PFBlockElement>& elements,
3443 const reco::PFBlock::LinkData& linkData,
3444 std::vector<bool>& active,
3445 std::vector<double>& psEne) {
3446
3447
3448
3449
3450
3451
3452 std::multimap<double, unsigned> sortedPS;
3453 block.associatedElements(iEcal, linkData, sortedPS, psElementType, reco::PFBlock::LINKTEST_ALL);
3454
3455
3456 for (auto const& ps : sortedPS) {
3457
3458 unsigned iPS = ps.second;
3459
3460
3461
3462 if (!active[iPS])
3463 continue;
3464
3465
3466 std::multimap<double, unsigned> sortedECAL;
3467 block.associatedElements(iPS, linkData, sortedECAL, reco::PFBlockElement::ECAL, reco::PFBlock::LINKTEST_ALL);
3468 unsigned jEcal = sortedECAL.begin()->second;
3469 if (jEcal != iEcal)
3470 continue;
3471
3472
3473 PFBlockElement::Type pstype = elements[iPS].type();
3474 assert(pstype == psElementType);
3475 PFClusterRef psclusterref = elements[iPS].clusterRef();
3476 assert(!psclusterref.isNull());
3477 psEne[0] += psclusterref->energy();
3478 active[iPS] = false;
3479 }
3480 }
3481
3482 bool PFAlgo::isFromSecInt(const reco::PFBlockElement& eTrack, string order) const {
3483 reco::PFBlockElement::TrackType T_TO_DISP = reco::PFBlockElement::T_TO_DISP;
3484 reco::PFBlockElement::TrackType T_FROM_DISP = reco::PFBlockElement::T_FROM_DISP;
3485
3486 reco::PFBlockElement::TrackType T_FROM_V0 = reco::PFBlockElement::T_FROM_V0;
3487
3488 bool bPrimary = (order.find("primary") != string::npos);
3489 bool bSecondary = (order.find("secondary") != string::npos);
3490 bool bAll = (order.find("all") != string::npos);
3491
3492 bool isToDisp = usePFNuclearInteractions_ && eTrack.trackType(T_TO_DISP);
3493 bool isFromDisp = usePFNuclearInteractions_ && eTrack.trackType(T_FROM_DISP);
3494
3495 if (bPrimary && isToDisp)
3496 return true;
3497 if (bSecondary && isFromDisp)
3498 return true;
3499 if (bAll && (isToDisp || isFromDisp))
3500 return true;
3501
3502
3503
3504
3505
3506 bool isFromDecay = (bAll || bSecondary) && usePFDecays_ && eTrack.trackType(T_FROM_V0);
3507
3508 return isFromDecay;
3509 }
3510
3511 void PFAlgo::postCleaning() {
3512
3513 double metX = 0.;
3514 double metY = 0.;
3515 double sumet = 0;
3516 std::vector<unsigned int> pfCandidatesToBeRemoved;
3517 for (auto const& pfc : *pfCandidates_) {
3518 metX += pfc.px();
3519 metY += pfc.py();
3520 sumet += pfc.pt();
3521 }
3522 double met2 = metX * metX + metY * metY;
3523
3524 double significance = std::sqrt(met2 / sumet);
3525 double significanceCor = significance;
3526 if (significance > minSignificance_) {
3527 double metXCor = metX;
3528 double metYCor = metY;
3529 double sumetCor = sumet;
3530 double met2Cor = met2;
3531 double deltaPhi = 3.14159;
3532 double deltaPhiPt = 100.;
3533 bool next = true;
3534 unsigned iCor = 1E9;
3535
3536
3537 while (next) {
3538 double metReduc = -1.;
3539
3540 for (unsigned i = 0; i < pfCandidates_->size(); ++i) {
3541 const PFCandidate& pfc = (*pfCandidates_)[i];
3542
3543
3544 if (pfc.particleId() != reco::PFCandidate::h_HF && pfc.particleId() != reco::PFCandidate::egamma_HF)
3545 continue;
3546
3547
3548 if (pfc.pt() < minHFCleaningPt_)
3549 continue;
3550
3551
3552 const bool skip = std::any_of(
3553 pfCandidatesToBeRemoved.begin(), pfCandidatesToBeRemoved.end(), [&](unsigned int j) { return i == j; });
3554 if (skip)
3555 continue;
3556
3557
3558 deltaPhi = std::acos((metX * pfc.px() + metY * pfc.py()) / (pfc.pt() * std::sqrt(met2)));
3559 deltaPhiPt = deltaPhi * pfc.pt();
3560 if (deltaPhiPt > maxDeltaPhiPt_)
3561 continue;
3562
3563
3564 double metXInt = metX - pfc.px();
3565 double metYInt = metY - pfc.py();
3566 double sumetInt = sumet - pfc.pt();
3567 double met2Int = metXInt * metXInt + metYInt * metYInt;
3568 if (met2Int < met2Cor) {
3569 metXCor = metXInt;
3570 metYCor = metYInt;
3571 metReduc = (met2 - met2Int) / met2Int;
3572 met2Cor = met2Int;
3573 sumetCor = sumetInt;
3574 significanceCor = std::sqrt(met2Cor / sumetCor);
3575 iCor = i;
3576 }
3577 }
3578
3579
3580 if (metReduc > minDeltaMet_) {
3581 pfCandidatesToBeRemoved.push_back(iCor);
3582 metX = metXCor;
3583 metY = metYCor;
3584 sumet = sumetCor;
3585 met2 = met2Cor;
3586 } else {
3587
3588 next = false;
3589 }
3590 }
3591
3592
3593 if (significance - significanceCor > minSignificanceReduction_ && significanceCor < maxSignificance_) {
3594 edm::LogInfo("PFAlgo|postCleaning") << "Significance reduction = " << significance << " -> " << significanceCor
3595 << " = " << significanceCor - significance;
3596 for (unsigned int toRemove : pfCandidatesToBeRemoved) {
3597 edm::LogInfo("PFAlgo|postCleaning") << "Removed : " << (*pfCandidates_)[toRemove];
3598 pfCleanedCandidates_.push_back((*pfCandidates_)[toRemove]);
3599 (*pfCandidates_)[toRemove].rescaleMomentum(1E-6);
3600
3601
3602 }
3603 }
3604 }
3605 }
3606
3607 void PFAlgo::checkCleaning(const reco::PFRecHitCollection& cleanedHits) {
3608
3609 if (cleanedHits.empty())
3610 return;
3611
3612
3613 double metX = 0.;
3614 double metY = 0.;
3615 double sumet = 0;
3616 std::vector<unsigned int> hitsToBeAdded;
3617 for (auto const& pfc : *pfCandidates_) {
3618 metX += pfc.px();
3619 metY += pfc.py();
3620 sumet += pfc.pt();
3621 }
3622 double met2 = metX * metX + metY * metY;
3623 double met2_Original = met2;
3624
3625
3626
3627 double metXCor = metX;
3628 double metYCor = metY;
3629 double sumetCor = sumet;
3630 double met2Cor = met2;
3631 bool next = true;
3632 unsigned iCor = 1E9;
3633
3634
3635 while (next) {
3636 double metReduc = -1.;
3637
3638 for (unsigned i = 0; i < cleanedHits.size(); ++i) {
3639 const PFRecHit& hit = cleanedHits[i];
3640 double length = std::sqrt(hit.position().mag2());
3641 double px = hit.energy() * hit.position().x() / length;
3642 double py = hit.energy() * hit.position().y() / length;
3643 double pt = std::sqrt(px * px + py * py);
3644
3645
3646 bool skip = false;
3647 for (unsigned int hitIdx : hitsToBeAdded) {
3648 if (i == hitIdx)
3649 skip = true;
3650 if (skip)
3651 break;
3652 }
3653 if (skip)
3654 continue;
3655
3656
3657 double metXInt = metX + px;
3658 double metYInt = metY + py;
3659 double sumetInt = sumet + pt;
3660 double met2Int = metXInt * metXInt + metYInt * metYInt;
3661
3662
3663 if (met2Int < met2Cor) {
3664 metXCor = metXInt;
3665 metYCor = metYInt;
3666 metReduc = (met2 - met2Int) / met2Int;
3667 met2Cor = met2Int;
3668 sumetCor = sumetInt;
3669
3670 iCor = i;
3671 }
3672 }
3673
3674
3675
3676 if (metReduc > minDeltaMet_) {
3677 hitsToBeAdded.push_back(iCor);
3678 metX = metXCor;
3679 metY = metYCor;
3680 sumet = sumetCor;
3681 met2 = met2Cor;
3682 } else {
3683
3684 next = false;
3685 }
3686 }
3687
3688
3689 if (std::sqrt(met2_Original) - std::sqrt(met2) > 5.) {
3690 LogTrace("PFAlgo|checkCleaning") << hitsToBeAdded.size() << " hits were re-added ";
3691 LogTrace("PFAlgo|checkCleaning") << "MET reduction = " << std::sqrt(met2_Original) << " -> " << std::sqrt(met2Cor)
3692 << " = " << std::sqrt(met2Cor) - std::sqrt(met2_Original);
3693 LogTrace("PFAlgo|checkCleaning") << "Added after cleaning check : ";
3694 for (unsigned int hitIdx : hitsToBeAdded) {
3695 const PFRecHit& hit = cleanedHits[hitIdx];
3696 PFCluster cluster(hit.layer(), hit.energy(), hit.position().x(), hit.position().y(), hit.position().z());
3697 reconstructCluster(cluster, hit.energy());
3698 LogTrace("PFAlgo|checkCleaning") << pfCandidates_->back() << ". time = " << hit.time();
3699 }
3700 }
3701 }