ATLAS Offline Software
Loading...
Searching...
No Matches
PFSubtractionTool.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2026 CERN for the benefit of the ATLAS collaboration
3*/
4
5#include "PFSubtractionTool.h"
6
7#include "eflowCaloObject.h"
11#include "eflowRecTrack.h"
14#include "PFClusterFiller.h"
15#include "PFTrackFiller.h"
16
20
21using namespace eflowSubtract;
22
23PFSubtractionTool::PFSubtractionTool(const std::string &type, const std::string &name, const IInterface *parent) : base_class(type, name, parent),
25{
26}
27
29= default;
30
32{
33
34 ATH_CHECK(m_theEOverPTool.retrieve());
35
36 ATH_CHECK(m_theEOverPTool->fillBinnedParameters(m_binnedParameters.get()));
37
39 if (!m_trkpos)
40 {
41 ATH_MSG_ERROR("Failed to get TrackPositionProvider for cluster preselection!");
42 return StatusCode::FAILURE;
43 }
44
45 //Retrieve track-cluster matching tools
46 ATH_CHECK(m_theMatchingTool.retrieve());
49
51
52 //Set the level of the helpers to the same as the tool here
53 m_pfSubtractionStatusSetter.msg().setLevel(this->msg().level());
54 m_pfSubtractionEnergyRatioCalculator.msg().setLevel(this->msg().level());
55 m_subtractor.m_facilitator.msg().setLevel(this->msg().level());
56
59 }
60
62
63 return StatusCode::SUCCESS;
64}
65
66void PFSubtractionTool::execute(const EventContext& ctx, eflowCaloObjectContainer *theEflowCaloObjectContainer, eflowRecTrackContainer *recTrackContainer, eflowRecClusterContainer *recClusterContainer) const
67{
68
69 ATH_MSG_DEBUG("Executing");
70
72 data.caloObjects = theEflowCaloObjectContainer;
73 const PFTrackFiller pfTrackFiller;
76
77 const PFClusterFiller pfClusterFiller;
80
81 ATH_MSG_DEBUG("This event has " << data.tracks.size() << " tracks " << data.clusters.size() << " clusters ");
82
83 unsigned int numMatches = matchAndCreateEflowCaloObj(ctx, data);
84
85 if (msgLvl(MSG::DEBUG)) printAllClusters(*recClusterContainer);
86
87 if (!m_calcEOverP){
90 else performSubtraction(numMatches,data);
91 }
93 }
94 //eoverp mode calculation has been moved to a dedicated PFBaseTool.
95
96}
97
98unsigned int PFSubtractionTool::matchAndCreateEflowCaloObj(const EventContext& ctx, PFData &data) const{
99
100 //Counts up how many tracks found at least 1 calorimeter cluster matched to it.
101 unsigned int nMatches(0);
102
103 /* Cache the original number of eflowCaloObjects, if there were any */
104 const unsigned int nCaloObj = data.caloObjects->size();
105
106 /* loop tracks in data.tracks and do matching */
107 for (auto *thisEfRecTrack : data.tracks)
108 {
110 if (!thisEfRecTrack->hasBin()) {
111 std::unique_ptr<eflowCaloObject> thisEflowCaloObject = std::make_unique<eflowCaloObject>();
112 thisEflowCaloObject->addTrack(thisEfRecTrack);
113 data.caloObjects->push_back(std::move(thisEflowCaloObject));
114 continue;
115 }
116
117 if (msgLvl(MSG::DEBUG))
118 {
119 const xAOD::TrackParticle *track = thisEfRecTrack->getTrack();
120 ATH_MSG_DEBUG("Matching track with e,pt, eta and phi " << track->e() << ", " << track->pt() << ", " << track->eta() << " and " << track->phi());
121 }
122
123 std::vector<eflowTrackClusterLink*> bestClusters;
124 std::vector<float> deltaRPrime;
125
127
128 const xAOD::TruthParticle* trackMatchedTruthParticle = nullptr;
130
131 const static SG::Accessor<TruthLink> truthLinkAccessor("truthParticleLink");
132
133 TruthLink truthLink = truthLinkAccessor(*(thisEfRecTrack->getTrack()));
134 //if not valid don't print a WARNING because this is an expected condition as discussed here:
135 //https://indico.cern.ch/event/795039/contributions/3391771/attachments/1857138/3050771/TruthTrackFTAGWS.pdf
136 if (truthLink.isValid()) trackMatchedTruthParticle = *truthLink;
137
138 if (trackMatchedTruthParticle){
139 double uniqueID = HepMC::uniqueID(trackMatchedTruthParticle);
140
142 if (!caloClusterReadDecorHandleNLeadingTruthParticles.isValid()){
143 ATH_MSG_WARNING("Failed to retrieve CaloCluster decoration with key " << caloClusterReadDecorHandleNLeadingTruthParticles.key());
144 }
145
146 for (auto * thisCluster : data.clusters){
147 //accessor for decoration
148 //split key into substring to get the name of the decoration
149
150 std::string decorHandleName = m_caloClusterReadDecorHandleKeyNLeadingTruthParticles.key();
151 std::string::size_type pos = decorHandleName.find(".");
152 std::string decorName = decorHandleName.substr(pos+1);
153
155
156 std::vector<std::pair<unsigned int, double > > uniqueIDTruthPairs = accessor(*(thisCluster->getCluster()));
157
158 for (auto &uniqueIDTruthPair : uniqueIDTruthPairs){
159 if (uniqueIDTruthPair.first == uniqueID){
160 eflowTrackClusterLink* thisLink = eflowTrackClusterLink::getInstance(thisEfRecTrack, thisCluster, ctx);
161 bestClusters.push_back(thisLink);
162 break;
163 }
164 }//loop over calocluster truth pair decorations
165 }//loop over caloclusters
166
167 }//if have truth particle matched to track
168 else ATH_MSG_VERBOSE("Track with pt, eta and phi " << thisEfRecTrack->getTrack()->pt() << ", " << thisEfRecTrack->getTrack()->eta() << " and " << thisEfRecTrack->getTrack()->phi() << " does not have a valid truth pointer");
169 }
170 else if (!m_recoverSplitShowers){
175 std::vector<std::pair<eflowRecCluster *, float>> bestClusters_02 = m_theMatchingToolForPull_02->doMatches(thisEfRecTrack, data.clusters, -1);
176 for (auto &matchpair : bestClusters_02)
177 {
178 eflowRecCluster *theCluster = matchpair.first;
179 float distancesq = matchpair.second;
180 eflowTrackClusterLink *trackClusterLink = eflowTrackClusterLink::getInstance(thisEfRecTrack, theCluster, ctx);
181 if (distancesq < 0.15 * 0.15)
182 {
183 // Narrower cone is a subset of the selected clusters
184 // Distance returned is deltaR^2
185 thisEfRecTrack->addAlternativeClusterMatch(trackClusterLink, "cone_015");
186 }
187 thisEfRecTrack->addAlternativeClusterMatch(trackClusterLink, "cone_02");
188 }//loop over bestClusters_02
189
190 //This matching scheme is used to match the calorimeter cluster(s) to be used in the charged showers subtraction for this track.
191 std::vector<std::pair<eflowRecCluster *, float>> matchedClusters = m_theMatchingTool->doMatches(thisEfRecTrack, data.clusters,m_nClusterMatchesToUse);
192 for (const auto& thePair : matchedClusters) {
193 bestClusters.push_back(eflowTrackClusterLink::getInstance(thisEfRecTrack, thePair.first, ctx));
194 if (m_addCPData) deltaRPrime.push_back(std::sqrt(thePair.second));
195 }
196 }
197 else {
198 const std::vector<eflowTrackClusterLink*>* matchedClusters_02 = thisEfRecTrack->getAlternativeClusterMatches("cone_02");
199 if (!matchedClusters_02) continue;
200 else bestClusters = *matchedClusters_02;
201 }
202
203 if (bestClusters.empty()) continue;
204
205 if (msgLvl(MSG::DEBUG))
206 {
207 for (auto *thisClusterLink : bestClusters ) {
208 xAOD::CaloCluster* thisCluster = thisClusterLink->getCluster()->getCluster();
209 ATH_MSG_DEBUG("Matched this track to cluster with e,pt, eta and phi " << thisCluster->e() << ", " << thisCluster->pt() << ", " << thisCluster->eta() << " and " << thisCluster->phi());
210 }
211 }
212
213 nMatches++;
214
215 //loop over the matched calorimeter clusters and associate tracks and clusters to each other as needed.
216 unsigned int linkIndex = 0;
217 for (auto *trkClusLink : bestClusters){
218
219 eflowRecCluster *thisEFRecCluster = trkClusLink->getCluster();
220
222 // Look up whether this cluster is intended for recovery
223 if (std::find(data.clusters.begin(), data.clusters.end(), trkClusLink->getCluster()) == data.clusters.end()) {
224 linkIndex++;
225 continue;
226 }
227 }
228
229 eflowTrackClusterLink *trackClusterLink = eflowTrackClusterLink::getInstance(thisEfRecTrack, thisEFRecCluster, ctx);
230 thisEfRecTrack->addClusterMatch(trackClusterLink);
232 thisEfRecTrack->addDeltaRPrime(deltaRPrime[linkIndex]);
233 }
234 thisEFRecCluster->addTrackMatch(trackClusterLink);
235 }
236 linkIndex++;
237 }
238
239 /* Create 3 types eflowCaloObjects: track-only, cluster-only, track-cluster-link */
240 std::vector<eflowRecCluster *> clusters(data.clusters.begin(), data.clusters.end());
241 if (m_recoverSplitShowers) std::sort(clusters.begin(), clusters.end(), eflowRecCluster::SortDescendingPt());
242 unsigned int nCaloObjects = eflowCaloObjectMaker::makeTrkCluCaloObjects(data.tracks, clusters, data.caloObjects);
243 ATH_MSG_DEBUG("Created " << nCaloObjects << " eflowCaloObjects.");
244 if (msgLvl(MSG::DEBUG)){
245 for (auto thisEFlowCaloObject : *(data.caloObjects)){
246 ATH_MSG_DEBUG("This eflowCaloObject has " << thisEFlowCaloObject->nTracks() << " tracks and " << thisEFlowCaloObject->nClusters() << " clusters ");
247 for (unsigned int count = 0; count < thisEFlowCaloObject->nTracks(); count++){
248 const xAOD::TrackParticle* thisTrack = thisEFlowCaloObject->efRecTrack(count)->getTrack();
249 ATH_MSG_DEBUG("Have track with e, pt, eta and phi of " << thisTrack->e() << ", " << thisTrack->pt() << ", " << thisTrack->eta() << " and " << thisTrack->phi());
250 }
251 for (unsigned int count = 0; count < thisEFlowCaloObject->nClusters(); count++){
252 const xAOD::CaloCluster* thisCluster = thisEFlowCaloObject->efRecCluster(count)->getCluster();
253 ATH_MSG_DEBUG("Have cluster with e, pt, eta and phi of " << thisCluster->e() << ", " << thisCluster->pt() << ", " << thisCluster->eta() << " and " << thisCluster->phi());
254 }
255 }
256 }
257
258 const double gaussianRadius = 0.032;
259 const double gaussianRadiusError = 1.0e-3;
260 const double maximumRadiusSigma = 3.0;
261
262 eflowLayerIntegrator integrator(gaussianRadius, gaussianRadiusError, maximumRadiusSigma, m_isHLLHC);
263
265 if (!m_recoverSplitShowers && 0 != nCaloObj) ATH_MSG_WARNING("Not in Split Showers Mode and already have " << nCaloObj << " eflowCaloObjects");
266
267 //For each eflowCaloObject we calculate the expected energy deposit in the calorimeter and cell ordering for subtraction.
268 for (unsigned int iCalo = nCaloObj; iCalo < data.caloObjects->size(); ++iCalo) {
269 eflowCaloObject* thisEflowCaloObject = data.caloObjects->at(iCalo);
270 thisEflowCaloObject->simulateShower(ctx, &integrator, m_binnedParameters.get(), m_useNNEnergy ? &(*m_NNEnergyPredictorTool) : nullptr, m_useLegacyEBinIndex);
271 if (m_useTruthForChargedShowerSubtraction) m_theTruthShowerSimulator->simulateShower(*thisEflowCaloObject);
272
273 }
274
275 if (!m_recoverSplitShowers) return nMatches;
276 else return nCaloObj;
277}
278
279void PFSubtractionTool::performSubtraction(unsigned int startingPoint,PFData &data) const{
280 unsigned int nEFCaloObs = data.caloObjects->size();
281 for (unsigned int iCalo = startingPoint; iCalo < nEFCaloObs; ++iCalo) {
282 eflowCaloObject* thisEflowCaloObject = data.caloObjects->at(iCalo);
283 this->performSubtraction(*thisEflowCaloObject);
284 }
285}
286
288
289 ATH_MSG_DEBUG("In performSubtraction");
290
291 unsigned int nClusters = thisEflowCaloObject.nClusters();
292 unsigned int nTrackMatches = thisEflowCaloObject.nTracks();
293
294 ATH_MSG_DEBUG("Have got an eflowCaloObject with " << nClusters << " clusters and " << nTrackMatches << " track matches");
295
296 if (msgLevel(MSG::DEBUG)){
297 for (unsigned int iTrack = 0; iTrack < nTrackMatches; ++iTrack){
298 eflowRecTrack* thisTrack = thisEflowCaloObject.efRecTrack(iTrack);
299 ATH_MSG_DEBUG("eflowCaloObject has track with E, pt and eta " << thisTrack->getTrack()->e() << ", " << thisTrack->getTrack()->pt() << " and " << thisTrack->getTrack()->eta());
300 }
301 }
302
303 //To keep historical behaviour when in recover split showers mode allow tracks with no cluster matches to proceed.
304 if (!m_recoverSplitShowers && nClusters < 1) return;
305
306 //Need at least one track in this eflowCaloObject to continue.
307 if (nTrackMatches < 1) return;
308
309 double expectedEnergy = thisEflowCaloObject.getExpectedEnergy();
310 double clusterEnergy = thisEflowCaloObject.getClusterEnergy();
311 double expectedSigma = std::sqrt(thisEflowCaloObject.getExpectedVariance());
312
313 /* Check e/p, if on first pass - return if e/p not consistent with expected e/p */
315 if (isEOverPFail(expectedEnergy, expectedSigma, clusterEnergy)) return;
316 }
317
318 const std::vector<std::pair<eflowTrackClusterLink *, std::pair<float, float>>> &matchedTrackList = thisEflowCaloObject.efRecLink();
319
320 ATH_MSG_DEBUG("Matched Track List has size " << matchedTrackList.size());
321
322 if (msgLevel(MSG::DEBUG))
323 {
324 for (unsigned int iTrack = 0; iTrack < nTrackMatches; ++iTrack)
325 {
326 const xAOD::TrackParticle *thisTrack = thisEflowCaloObject.efRecTrack(iTrack)->getTrack();
327 ATH_MSG_DEBUG("eflowCaloObject has track match with E, pt and eta " << thisTrack->e() << ", " << thisTrack->pt() << " and " << thisTrack->eta());
328 }
329 }
330
331 ATH_MSG_DEBUG("About to perform subtraction for this eflowCaloObject");
332
333 bool wasAnnihilated = false;
334
335 //First deal with non-split showers mode
337 /* Check if we can annihilate right away - true if matched cluster has only the expected energy deposit */
338 if (canAnnihilate(expectedEnergy, expectedSigma, clusterEnergy)){
339
340 wasAnnihilated = true;
341
342 std::vector<std::pair<xAOD::CaloCluster *, bool>> clusterList;
343 std::map<xAOD::CaloCluster *, double> clusterEnergyMap;
344 unsigned nCluster = thisEflowCaloObject.nClusters();
345 for (unsigned iCluster = 0; iCluster < nCluster; ++iCluster){
346 clusterList.emplace_back(thisEflowCaloObject.efRecCluster(iCluster)->getCluster(), false);
347 }
348
349 ATH_MSG_DEBUG("We are going to annihilate. ExpectedEnergy, expectedSigma and clusterEnergy are " << expectedEnergy << ", " << expectedSigma << " and " << clusterEnergy);
350 if (msgLevel(MSG::DEBUG))
351 for (const auto& thisPair : clusterList)
352 ATH_MSG_DEBUG("Annihilating cluster with E and eta " << thisPair.first->e() << " and " << thisPair.first->eta());
353
354 m_pfSubtractionStatusSetter.markAllTracksAnnihStatus(thisEflowCaloObject);
355
356 //before we remove all the cells, we create a list of the removed cells if in doCPData mode
357 if (m_addCPData) this->addSubtractedCells(thisEflowCaloObject, clusterList);
358
360
361 if (msgLevel(MSG::DEBUG))
362 for (const auto& thisPair : clusterList)
363 ATH_MSG_DEBUG("Have Annihilated cluster with E and eta " << thisPair.first->e() << " and " << thisPair.first->eta());
364
365 /* Flag all tracks in this system as subtracted */
366 for (unsigned iTrack = 0; iTrack < thisEflowCaloObject.nTracks(); ++iTrack){
367 eflowRecTrack *thisEfRecTrack = (matchedTrackList[iTrack].first)->getTrack();
368 if (!thisEfRecTrack->isSubtracted()) thisEfRecTrack->setSubtracted();
369 }
370
371 }//if can annihilate this track-cluster systems matched cluster
372 }//split shower recovery mode or regular mode where above annihilation was not triggered
373 if (m_recoverSplitShowers || !wasAnnihilated){
374
375 for (unsigned iTrack = 0; iTrack < thisEflowCaloObject.nTracks(); ++iTrack){
376
377 eflowRecTrack *thisEfRecTrack = thisEflowCaloObject.efRecTrack(iTrack);
378
379 ATH_MSG_DEBUG("About to subtract track with e, pt, eta and phi of " << thisEfRecTrack->getTrack()->e() << ", " << thisEfRecTrack->getTrack()->pt() << ", " << thisEfRecTrack->getTrack()->eta() << " and "
380 << thisEfRecTrack->getTrack()->eta());
381
382 if (!thisEfRecTrack->hasBin()) continue;
383
384 ATH_MSG_DEBUG("Have bin for this eflowCaloObject");
385
386 if (thisEfRecTrack->isInDenseEnvironment() && !m_recoverSplitShowers) continue;
387
388 ATH_MSG_DEBUG("Am not in dense environment for this eflowCaloObject");
389
390 /* Get matched cluster via Links */
391 std::vector<eflowRecCluster *> matchedClusters;
392 const std::vector<eflowTrackClusterLink *>& links = thisEfRecTrack->getClusterMatches();
393 matchedClusters.reserve(links.size());
394 for (auto* thisEFlowTrackClusterLink : links)
395 matchedClusters.push_back(thisEFlowTrackClusterLink->getCluster());
397 std::sort(matchedClusters.begin(),
398 matchedClusters.end(),
400
401 if (msgLvl(MSG::DEBUG)) {
402 for (auto* thisClus : matchedClusters)
404 "Haved matched cluster "
405 << thisClus->getCluster()->index() << " with e,pt, eta and phi of "
406 << thisClus->getCluster()->e() << ", "
407 << thisClus->getCluster()->pt() << ", "
408 << thisClus->getCluster()->eta() << " and "
409 << thisClus->getCluster()->phi() << " will be subtracted");
410 }
411
412 /* Do subtraction */
413 std::vector<std::pair<xAOD::CaloCluster *, bool>> clusterSubtractionList;
414 clusterSubtractionList.reserve(matchedClusters.size());
415 std::map<xAOD::CaloCluster *, double> clusterEnergyMap;
416 for (auto *thisEFlowRecCluster : matchedClusters){
417 xAOD::CaloCluster *thisCluster = thisEFlowRecCluster->getCluster();
418 clusterSubtractionList.emplace_back(thisCluster, false);
419 clusterEnergyMap[thisCluster] = thisCluster->e();
420 }
421
422 ATH_MSG_DEBUG("Have filled clusterSubtractionList for this eflowCaloObject");
423
424 unsigned int trackIndex = thisEfRecTrack->getTrack()->index();
425
426 //Previously we only checked this in recover split showers, but makes sense to check it in both passes.
427 auto sumClusEnergy = [](double accumulator, std::pair<xAOD::CaloCluster *, bool> thisPair){ return accumulator += thisPair.first->e();};
428 double totalClusterEnergy = std::accumulate(clusterSubtractionList.begin(),clusterSubtractionList.end(),0.0,sumClusEnergy);
429
430 /* Check if we can annihilate right away - true if matched cluster has only the expected energy deposit */
431 if(canAnnihilate(thisEfRecTrack->getEExpect(),std::sqrt(thisEfRecTrack->getVarEExpect()),totalClusterEnergy)){
432
433 if (msgLevel(MSG::DEBUG))
434 for (const auto& thisPair : clusterSubtractionList)
435 ATH_MSG_DEBUG("Annihilating cluster with E and eta " << thisPair.first->e() << " and " << thisPair.first->eta());
436
437 //before we remove all the cells, we create a list of the removed cells if in doCPData mode
438 if (m_addCPData) this->addSubtractedCells(thisEflowCaloObject, clusterSubtractionList);
439
440 Subtractor::annihilateClusters(clusterSubtractionList);
441 //Now we should mark all of these clusters as being subtracted
442 //Now need to mark which clusters were modified in the subtraction procedure
443 std::vector<std::pair<float, float>> clusterSubtractedEnergyRatios;
444 m_pfSubtractionEnergyRatioCalculator.calculateSubtractedEnergyRatiosForAnnih(clusterSubtractionList, clusterEnergyMap, clusterSubtractedEnergyRatios);
445 m_pfSubtractionStatusSetter.markSubtractionStatus(clusterSubtractionList, clusterSubtractedEnergyRatios, thisEflowCaloObject, trackIndex);
446 }
447 else
448 {
449
450 /* Subtract the track from all matched clusters */
451 m_subtractor.subtractTracksFromClusters(thisEfRecTrack, clusterSubtractionList, m_addCPData);
452
453 //recalculate total cluster energy from the clusters afer subtraction
454 totalClusterEnergy = std::accumulate(clusterSubtractionList.begin(),clusterSubtractionList.end(),0.0,sumClusEnergy);
455
456 /* Annihilate the cluster(s) if the remnant is small (i.e. below k*sigma) */
457 if (canAnnihilate(0.0,std::sqrt(thisEfRecTrack->getVarEExpect()), totalClusterEnergy)){
458
459 if (msgLevel(MSG::DEBUG))
460 for (const auto& thisPair : clusterSubtractionList){
461 ATH_MSG_DEBUG("Annihilating remnant cluster with E and eta " << thisPair.first->e() << " and " << thisPair.first->eta());
462 }
463 eflowSubtract::Subtractor::annihilateClusters(clusterSubtractionList);
464 //Now we should mark all of these clusters as being subtracted
465 std::vector<std::pair<float, float>> clusterSubtractedEnergyRatios;
466 m_pfSubtractionEnergyRatioCalculator.calculateSubtractedEnergyRatiosForAnnih(clusterSubtractionList, clusterEnergyMap, clusterSubtractedEnergyRatios);
467 m_pfSubtractionStatusSetter.markSubtractionStatus(clusterSubtractionList, clusterSubtractedEnergyRatios, thisEflowCaloObject, trackIndex);
468 }//if remove the remnant after cell by cell subtraction
469 else
470 {
471 std::vector<std::pair<float, float>> clusterSubtractedEnergyRatios;
472 m_pfSubtractionEnergyRatioCalculator.calculateSubtractedEnergyRatios(clusterSubtractionList, clusterEnergyMap, clusterSubtractedEnergyRatios);
473 m_pfSubtractionStatusSetter.markSubtractionStatus(clusterSubtractionList, clusterSubtractedEnergyRatios, thisEflowCaloObject, trackIndex);
474 }//if don't remove the remnant after cell by cell subtraction
475
476 }//if not annihilating, and instead subtracting cell by cell
477
478 ATH_MSG_DEBUG("Have subtracted charged shower for this eflowRecTrack");
479
480 /* Flag tracks as subtracted */
481 if (!thisEfRecTrack->isSubtracted()) thisEfRecTrack->setSubtracted();
482
483 }//loop over tracks in eflowCaloObject
484 }//cell by cell subtraction
485
486}
487
489
490 ATH_MSG_DEBUG("In performTruthSubtraction");
491
492 unsigned int nEFCaloObs = data.caloObjects->size();
493
494 for (unsigned int iCalo = 0; iCalo < nEFCaloObs; ++iCalo) {
495 eflowCaloObject* thisEflowCaloObject = data.caloObjects->at(iCalo);
496 this->performTruthSubtraction(*thisEflowCaloObject);
497 }
498
499}
500
502
503 for (unsigned iTrack = 0; iTrack < thisEflowCaloObject.nTracks(); ++iTrack){
504 eflowRecTrack *thisEfRecTrack = thisEflowCaloObject.efRecTrack(iTrack);
505
506 //although we are subtracting the truth, to be consistent we only do it if a reco
507 //e/p lookup bin exists for this track
508 if (!thisEfRecTrack->hasBin()) continue;
509
510 //Similarly we skip tracks in a dense environment
511 if (thisEfRecTrack->isInDenseEnvironment()) continue;
512
513 thisEfRecTrack->setSubtracted();
514
515 //get the set of matched clusters
516 const std::vector<eflowTrackClusterLink *>& links = thisEfRecTrack->getClusterMatches();
517
518 for (auto thisLink : links){
519 xAOD::CaloCluster *thisCluster = thisLink->getCluster()->getCluster();
520 CaloClusterCellLink* theCellLinks = thisCluster->getOwnCellLinks();
521 CaloClusterCellLink::iterator theCell = theCellLinks->begin();
522 CaloClusterCellLink::iterator lastCell = theCellLinks->end();
523
524 //loop over the cells in this cluster and subtract shower using truth information
525 //We can either remove a cell entireley if it has any truth deposit (closer to what the real
526 //reco algorithm does) or reweight the cells contribution based on subtracting the truth
527 //energy from the reco cell energy
528 //We only advance the iterator, theCell, if we *dont* call removeCell to avoid issues with
529 //invalid iterators
530 //We also have to reset the lastCell iterator after each call to ensure the loop exits at the end,
531 //instead of being stuck in an infinite loop.
532 for (; theCell != lastCell;){
533 //get the truth energy for this cell
534 double truthEnergy = thisEfRecTrack->getCellTruthEnergy(*theCell);
535 //reweight the cell such that energy*weight gives the new energy
536 double oldCellEnergy = theCell->energy()*(theCell.weight());
537 double subtractedCellWeight = (oldCellEnergy - truthEnergy)/oldCellEnergy;
538
539 if (0.0 != truthEnergy && m_useFullCellTruthSubtraction) {
540 thisCluster->removeCell(*theCell);
541 lastCell = theCellLinks->end();
542 }
544 theCell.reweight(subtractedCellWeight);
545 ++theCell;
546 }
547 else ++theCell;
548
549 }//cell loop
550
551 float oldEnergy = thisCluster->e();
552 CaloClusterKineHelper::calculateKine(thisCluster, true, true);
553 if (0.0 != oldEnergy) {
554 float energyAdjustment = thisCluster->e() / oldEnergy;
555 thisCluster->setRawE(thisCluster->rawE() * energyAdjustment);
556 thisCluster->setRawEta(thisCluster->eta());
557 thisCluster->setRawPhi(thisCluster->phi());
558 }
559 }
560
561 }//eflowCaloObject track loop
562
563}
564
565bool PFSubtractionTool::isEOverPFail(double expectedEnergy, double sigma, double clusterEnergy) const
566{
567 if ((expectedEnergy == 0) && (clusterEnergy > 0)) return false;
568 return clusterEnergy < expectedEnergy - m_consistencySigmaCut * sigma;
569}
570
571bool PFSubtractionTool::canAnnihilate(double expectedEnergy, double sigma, double clusterEnergy) const
572{
573 return clusterEnergy - expectedEnergy < m_subtractionSigmaCut * sigma;
574}
575
577 std::stringstream result;
578 result << " track with E, eta and phi "<< track->e() << ", " << track->eta() << " and " << track->phi();
579 return result.str();
580}
581
583 std::stringstream result;
584 result << " cluster with E, eta and phi of " << cluster->e() << ", " << cluster->eta() << " and " << cluster->phi();
585 return result.str();
586}
587
588void PFSubtractionTool::printAllClusters(const eflowRecClusterContainer& recClusterContainer) const {
589
590 for ( const auto *thisEFRecCluster : recClusterContainer){
591 if (thisEFRecCluster->getTrackMatches().empty()) {
592 ATH_MSG_DEBUG("Isolated" << printCluster(thisEFRecCluster->getCluster()));
593 } else {
594 ATH_MSG_DEBUG("Matched" << printCluster(thisEFRecCluster->getCluster()));
595 std::vector<eflowTrackClusterLink*> theTrackLinks = thisEFRecCluster->getTrackMatches();
596 for ( auto *thisTrack : theTrackLinks){
597 ATH_MSG_DEBUG("Matched" << printTrack(thisTrack->getTrack()->getTrack()));
598 }
599 }
600 }
601}
602
603void PFSubtractionTool::addSubtractedCells(eflowCaloObject& thisEflowCaloObject, const std::vector<std::pair<xAOD::CaloCluster *, bool> >& clusterList) const{
604
605 unsigned int numTracks = thisEflowCaloObject.nTracks();
606
607 for (unsigned int iTrack = 0; iTrack < numTracks; ++iTrack){
608 eflowRecTrack* thisTrack = thisEflowCaloObject.efRecTrack(iTrack);
609 for (const auto& thisPair : clusterList){
610 xAOD::CaloCluster* thisCluster = thisPair.first;
611 const CaloClusterCellLink* theCellLink = thisCluster->getCellLinks();
612 CaloClusterCellLink::const_iterator theCell = theCellLink->begin();
613 CaloClusterCellLink::const_iterator lastCell = theCellLink->end();
614 for (; theCell != lastCell; ++theCell) thisTrack->addSubtractedCaloCell(ElementLink<CaloCellContainer>("AllCalo",theCell.index()),theCell.weight()/numTracks);
615 }
616 }
617}
618
619StatusCode PFSubtractionTool::finalize() { return StatusCode::SUCCESS; }
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
#define ATH_MSG_DEBUG(x)
char data[hepevt_bytes_allocation_ATLAS]
Definition HepEvt.cxx:11
Handle class for reading a decoration on an object.
ElementLink< xAOD::TruthParticleContainer > TruthLink
double energy() const
get energy (data member)
Definition CaloCell.h:327
static void calculateKine(xAOD::CaloCluster *clu, const bool useweight=true, const bool updateLayers=true, const bool useGPUCriteria=false)
Helper class to calculate cluster kinematics based on cells.
static void fillClustersToRecover(PFData &data)
static void fillClustersToConsider(PFData &data, eflowRecClusterContainer &recClusterContainer)
static std::unique_ptr< IPositionProvider > Get(const std::string &positionType)
std::unique_ptr< PFMatch::TrackEtaPhiInFixedLayersProvider > m_trkpos
Track position provider to be used to preselect clusters.
Gaudi::Property< bool > m_useNNEnergy
Toggle whether we use the neural net energy.
PFSubtractionEnergyRatioCalculator m_pfSubtractionEnergyRatioCalculator
ToolHandle< IEFlowCellEOverPTool > m_theEOverPTool
Tool for getting e/p values and hadronic shower cell ordering principle parameters.
ToolHandle< PFTrackClusterMatchingTool > m_theMatchingTool
Default track-cluster matching tool.
void performTruthSubtraction(PFData &data) const
ToolHandle< PFEnergyPredictorTool > m_NNEnergyPredictorTool
Tool for getting predictiing the energy using an ONNX model.
Gaudi::Property< double > m_subtractionSigmaCut
Parameter that controls whether to use retain remaining calorimeter energy in track-cluster system,...
Gaudi::Property< double > m_consistencySigmaCut
Parameter that controls whether a track, in a track-cluster system, will be processed by the split sh...
void printAllClusters(const eflowRecClusterContainer &recClusterContainer) const
Gaudi::Property< bool > m_useFullCellTruthSubtraction
Toggle whether we fully remove a cell with a truth deposit or reweight it based on truth contribution...
void addSubtractedCells(eflowCaloObject &thisEflowCaloObject, const std::vector< std::pair< xAOD::CaloCluster *, bool > > &clusterList) const
eflowSubtract::Subtractor m_subtractor
SG::ReadDecorHandleKey< xAOD::CaloClusterContainer > m_caloClusterReadDecorHandleKeyNLeadingTruthParticles
Read handle key to decorate CaloCluster with threeN leading truth particle uniqueID and energy.
Gaudi::Property< int > m_nClusterMatchesToUse
Number of clusters to match to each track if not doing recover split shower subtraction.
PFSubtractionStatusSetter m_pfSubtractionStatusSetter
Gaudi::Property< bool > m_useTruthMatching
Toggle whether to cheat and use truth information for track-cluster matching - only for performance s...
bool canAnnihilate(double expectedEnergy, double sigma, double clusterEnergy) const
Gaudi::Property< bool > m_addCPData
Toggle whether to decorate eflowRecTrack with addutional data for Combined Performance studies.
void execute(const EventContext &ctx, eflowCaloObjectContainer *theEflowCaloObjectContainer, eflowRecTrackContainer *recTrackContainer, eflowRecClusterContainer *recClusterContainer) const
ToolHandle< PFSimulateTruthShowerTool > m_theTruthShowerSimulator
std::unique_ptr< eflowEEtaBinnedParameters > m_binnedParameters
PFSubtractionTool(const std::string &type, const std::string &name, const IInterface *parent)
static std::string printCluster(const xAOD::CaloCluster *cluster)
unsigned int matchAndCreateEflowCaloObj(const EventContext &ctx, PFData &data) const
This matches ID tracks and CaloClusters, and then creates eflowCaloObjects.
Gaudi::Property< bool > m_calcEOverP
Toggle EOverP algorithm mode, whereby no charged shower subtraction is performed.
Gaudi::Property< bool > m_useTruthForChargedShowerSubtraction
Toggle whether we use truth information for the charged shower subtraction or not.
ToolHandle< PFTrackClusterMatchingTool > m_theMatchingToolForPull_02
Gaudi::Property< bool > m_useLegacyEBinIndex
Further discussion about why this flag exists can be found in https://its.cern.ch/jira/browse/ATLJETM...
Gaudi::Property< bool > m_recoverSplitShowers
Toggle whether we are recovering split showers or not.
ToolHandle< PFTrackClusterMatchingTool > m_theMatchingToolForPull_015
static std::string printTrack(const xAOD::TrackParticle *track)
Gaudi::Property< bool > m_isHLLHC
Toggle whether we have the HLLHC setup.
void performSubtraction(unsigned int startingPoint, PFData &data) const
bool isEOverPFail(double expectedEnergy, double sigma, double clusterEnergy) const
static void fillTracksToConsider(PFData &data, eflowRecTrackContainer &recTrackContainer)
static void fillTracksToRecover(PFData &data)
Helper class to provide type-safe access to aux data.
Handle class for reading a decoration on an object.
static unsigned int makeTrkCluCaloObjects(eflowRecTrackContainer *eflowTrackContainer, eflowRecClusterContainer *eflowClusterContainer, eflowCaloObjectContainer *caloObjectContainer)
An internal EDM object which stores information about systems of associated tracks and calorimeter cl...
void simulateShower(const EventContext &ctx, eflowLayerIntegrator *integrator, const eflowEEtaBinnedParameters *binnedParameters, const PFEnergyPredictorTool *energyP, bool useLegacyEnergyBinIndexing)
double getExpectedVariance() const
double getClusterEnergy() const
unsigned nClusters() const
const eflowRecTrack * efRecTrack(int i) const
const std::vector< std::pair< eflowTrackClusterLink *, std::pair< float, float > > > & efRecLink() const
unsigned nTracks() const
double getExpectedEnergy() const
const eflowRecCluster * efRecCluster(int i) const
Inherits from eflowEEtaBinBase.
This class calculates the LHED (Layer of Highest Energy Density) in a cluster or group of clusters.
This class extends the information about a xAOD::CaloCluster.
void addTrackMatch(eflowTrackClusterLink *trackMatch)
xAOD::CaloCluster * getCluster()
This class extends the information about a xAOD::Track.
double getEExpect() const
const std::vector< eflowTrackClusterLink * > & getClusterMatches() const
bool isSubtracted() const
const xAOD::TrackParticle * getTrack() const
void addSubtractedCaloCell(ElementLink< CaloCellContainer > theCellLink, const double &weight)
bool hasBin() const
double getCellTruthEnergy(const CaloCell *cell) const
bool isInDenseEnvironment() const
double getVarEExpect() const
static void annihilateClusters(std::vector< std::pair< xAOD::CaloCluster *, bool > > &clusters)
void setRawEta(flt_t)
Set for signal state UNCALIBRATED.
flt_t rawE() const
void setRawPhi(flt_t)
Set for signal state UNCALIBRATED.
const CaloClusterCellLink * getCellLinks() const
Get a pointer to the CaloClusterCellLink object (const version).
virtual double pt() const
The transverse momentum ( ) of the particle (negative for negative-energy clusters).
void setRawE(flt_t)
Set Energy for signal state UNCALIBRATED.
virtual double eta() const
The pseudorapidity ( ) of the particle.
virtual double e() const
The total energy of the particle.
virtual double phi() const
The azimuthal angle ( ) of the particle.
CaloClusterCellLink * getOwnCellLinks()
Get a pointer to the owned CaloClusterCellLink object (non-const version).
bool removeCell(const CaloCell *ptr)
Method to remove a cell to the cluster (slow!) (Beware: Kinematics not updated!).
virtual double phi() const override final
The azimuthal angle ( ) of the particle (has range to .).
virtual double pt() const override final
The transverse momentum ( ) of the particle.
virtual double eta() const override final
The pseudorapidity ( ) of the particle.
virtual double e() const override final
The total energy of the particle.
static std::string release
Definition computils.h:50
int count(std::string s, const std::string &regx)
count how many occurances of a regx are in a string
Definition hcg.cxx:148
int uniqueID(const T &p)
STL namespace.
void sort(typename DataModel_detail::iterator< DVL > beg, typename DataModel_detail::iterator< DVL > end)
Specialization of sort for DataVector/List.
CaloCluster_v1 CaloCluster
Define the latest version of the calorimeter cluster class.
TrackParticle_v1 TrackParticle
Reference the current persistent version:
TruthParticle_v1 TruthParticle
Typedef to implementation.
MsgStream & msg
Definition testRead.cxx:32