9#include <ecl/modules/eclTrackClusterMatching/ECLTrackClusterMatchingModule.h>
11#include <ecl/dbobjects/ECLTrackClusterMatchingParameterizations.h>
12#include <ecl/dbobjects/ECLTrackClusterMatchingThresholds.h>
13#include <ecl/utility/utilityFunctions.h>
15#include <framework/gearbox/Const.h>
16#include <framework/logging/Logger.h>
17#include <mdst/dataobjects/HitPatternCDC.h>
25 m_matchingParameterizations(
"ECLTrackClusterMatchingParameterizations"),
26 m_matchingThresholds(
"ECLTrackClusterMatchingThresholds")
28 setDescription(
"Creates and saves a Relation between Tracks and ECLCluster in the DataStore. It uses the existing Relation between Tracks and ExtHit as well as the Relation between ECLCluster and ExtHit.");
31 "if true use track cluster matching based on angular distance, if false use matching based on entered crystals",
bool(
true));
33 "set false if you want to set the matching criterion on your own",
bool(
true));
35 "the 2D consistency of Delta theta and Delta phi has to exceed this value for a track to be matched to an ECL cluster", 1e-6);
37 "tracks with pt greater than this value will exclusively be matched based on angular distance", 0.3);
39 "distance of polar angle from gaps where crystal-entering based matching is applied (in rad)", 0.1);
41 "bad VXD-standalone tracks cause (too) low photon efficiency in end caps, temporarily fixed by requiring minimal number of CDC hits",
44 "switch to exclude tracks with zero charge from track-cluster matching",
bool(
false));
68 auto updateParameterizationFunctions = [
this]() {
92 auto updateMatchingThresholds = [
this]() {
99 updateParameterizationFunctions();
100 updateMatchingThresholds();
113 eclCluster.setIsTrack(
false);
125 set<int> uniqueShowerIds;
128 set<int> uniquehypothesisIds;
129 vector<int> hypothesisIds;
130 vector<double> energies;
131 vector<int> arrayIndexes;
137 for (
const auto& extHit : track.getRelationsTo<
ExtHit>()) {
139 const int cell = extHit.getCopyID() + 1;
143 [&](
const ECLCalDigit & d) { return d.getCellId() == cell; }
149 for (
auto& shower : idigit->getRelationsFrom<
ECLShower>()) {
150 bool inserted = (uniqueShowerIds.insert(shower.getUniqueId())).second;
153 if (!inserted)
continue;
155 hypothesisIds.push_back(shower.getHypothesisId());
156 energies.push_back(shower.getEnergy());
157 arrayIndexes.push_back(shower.getArrayIndex());
158 uniquehypothesisIds.insert(shower.getHypothesisId());
160 B2DEBUG(29, shower.getArrayIndex() <<
" " << shower.getHypothesisId() <<
" " << shower.getEnergy() <<
" " <<
161 shower.getConnectedRegionId());
167 for (
auto hypothesisId : uniquehypothesisIds) {
168 double highestEnergy = 0.0;
171 for (
unsigned ix = 0; ix < energies.size(); ix++) {
172 if (hypothesisIds[ix] == hypothesisId and energies[ix] > highestEnergy) {
173 highestEnergy = energies[ix];
174 arrayindex = arrayIndexes[ix];
179 if (arrayindex > -1) {
181 shower->setIsTrack(
true);
183 track.addRelationTo(shower);
184 track.addRelationTo(shower, 1.0,
"AngularDistance");
186 track.addRelationTo(shower, 1.0,
"EnterCrystal");
188 B2DEBUG(29, shower->getArrayIndex() <<
" " << shower->getIsTrack());
192 if (cluster !=
nullptr) {
193 cluster->setIsTrack(
true);
195 track.addRelationTo(cluster);
196 track.addRelationTo(cluster, 1.0,
"AngularDistance");
198 track.addRelationTo(cluster, 1.0,
"EnterCrystal");
206 if (track.getRelationsTo<
ECLCluster>(
"",
"AngularDistance").size() > 0)
continue;
210 ECL::DetectorRegion trackDetectorRegion = ECL::getDetectorRegion(theta);
213 map<int, pair<double, int>> hypothesisIdBestQualityCROSSArrayIndexMap;
214 map<int, pair<double, int>> hypothesisIdBestQualityDLArrayIndexMap;
215 map<int, pair<double, int>> hypothesisIdBestQualityNEARArrayIndexMap;
216 set<int> uniqueHypothesisIds;
221 for (
const auto& extHit : track.getRelationsTo<
ExtHit>()) {
224 if (!eclCluster)
continue;
226 if (eclShower !=
nullptr) {
229 if (abs(eclDetectorRegion - trackDetectorRegion) == 1)
continue;
232 double phiHit = extHit.getPosition().Phi();
233 double phiShower = eclShower->
getPhi();
234 double deltaPhi = phiHit - phiShower;
235 if (deltaPhi > M_PI) {
236 deltaPhi = deltaPhi - 2 * M_PI;
237 }
else if (deltaPhi < -M_PI) {
238 deltaPhi = deltaPhi + 2 * M_PI;
240 double thetaHit = extHit.getPosition().Theta();
241 double thetaShower = eclShower->
getTheta();
242 double deltaTheta = thetaHit - thetaShower;
244 double quality =
showerQuality(deltaPhi, deltaTheta, pt, eclDetectorRegion, extHitStatus);
246 bool inserted = (uniqueHypothesisIds.insert(hypothesisId)).second;
248 hypothesisIdBestQualityCROSSArrayIndexMap.insert(make_pair(hypothesisId, make_pair(0, -1)));
249 hypothesisIdBestQualityDLArrayIndexMap.insert(make_pair(hypothesisId, make_pair(0, -1)));
250 hypothesisIdBestQualityNEARArrayIndexMap.insert(make_pair(hypothesisId, make_pair(0, -1)));
252 if (extHitStatus == EXT_ECLCROSS) {
253 if (quality > hypothesisIdBestQualityCROSSArrayIndexMap.at(hypothesisId).first) {
254 hypothesisIdBestQualityCROSSArrayIndexMap[hypothesisId] = make_pair(quality, eclShower->
getArrayIndex());
256 }
else if (extHitStatus == EXT_ECLDL) {
257 if (quality > hypothesisIdBestQualityDLArrayIndexMap.at(hypothesisId).first) {
258 hypothesisIdBestQualityDLArrayIndexMap[hypothesisId] = make_pair(quality, eclShower->
getArrayIndex());
261 if (quality > hypothesisIdBestQualityNEARArrayIndexMap.at(hypothesisId).first) {
262 hypothesisIdBestQualityNEARArrayIndexMap[hypothesisId] = make_pair(quality, eclShower->
getArrayIndex());
268 vector<map<int, pair<double, int>>> hypothesisIdBestQualityArrayIndexMaps;
269 hypothesisIdBestQualityArrayIndexMaps.push_back(hypothesisIdBestQualityCROSSArrayIndexMap);
270 hypothesisIdBestQualityArrayIndexMaps.push_back(hypothesisIdBestQualityDLArrayIndexMap);
271 hypothesisIdBestQualityArrayIndexMaps.push_back(hypothesisIdBestQualityNEARArrayIndexMap);
273 for (
const auto& uniqueHypothesisId : uniqueHypothesisIds) {
274 for (
const auto& hypothesisIdBestQualityArrayIndexMap : hypothesisIdBestQualityArrayIndexMaps) {
276 && hypothesisIdBestQualityArrayIndexMap.at(uniqueHypothesisId).second > -1) {
277 auto shower =
m_eclShowers[hypothesisIdBestQualityArrayIndexMap.at(uniqueHypothesisId).second];
278 shower->setIsTrack(
true);
279 track.addRelationTo(shower);
280 track.addRelationTo(shower, 1.0,
"AngularDistance");
282 if (cluster !=
nullptr) {
283 cluster->setIsTrack(
true);
284 track.addRelationTo(cluster);
285 track.addRelationTo(cluster, 1.0,
"AngularDistance");
301 if ((extHit.
getDetectorID() != Const::EDetector::ECL))
return false;
302 if ((extHit.
getStatus() != EXT_ENTER))
return false;
303 if (extHit.
getCopyID() == -1)
return false;
309 if ((extHit.
getDetectorID() != Const::EDetector::ECL))
return false;
311 if (extHitStatus == EXT_ECLCROSS || extHitStatus == EXT_ECLDL || extHitStatus == EXT_ECLNEAR)
return true;
316 int eclDetectorRegion,
int hitStatus)
const
318 double phi_consistency =
phiConsistency(deltaPhi, pt, eclDetectorRegion, hitStatus);
319 double theta_consistency =
thetaConsistency(deltaTheta, pt, eclDetectorRegion, hitStatus);
320 return phi_consistency * theta_consistency * (1 - log(phi_consistency * theta_consistency));
326 if (eclDetectorRegion == ECL::DetectorRegion::FWD || eclDetectorRegion == ECL::DetectorRegion::FWDGap) {
327 if (hitStatus == EXT_ECLCROSS) {
329 }
else if (hitStatus == EXT_ECLDL) {
334 }
else if (eclDetectorRegion == ECL::DetectorRegion::BRL) {
335 if (hitStatus == EXT_ECLCROSS) {
337 }
else if (hitStatus == EXT_ECLDL) {
342 }
else if (eclDetectorRegion == ECL::DetectorRegion::BWD || eclDetectorRegion == ECL::DetectorRegion::BWDGap) {
343 if (hitStatus == EXT_ECLCROSS) {
345 }
else if (hitStatus == EXT_ECLDL) {
353 return erfc(abs(deltaPhi) / phi_RMS);
359 if (eclDetectorRegion == ECL::DetectorRegion::FWD || eclDetectorRegion == ECL::DetectorRegion::FWDGap) {
360 if (hitStatus == EXT_ECLCROSS) {
362 }
else if (hitStatus == EXT_ECLDL) {
367 }
else if (eclDetectorRegion == ECL::DetectorRegion::BRL) {
368 if (hitStatus == EXT_ECLCROSS) {
370 }
else if (hitStatus == EXT_ECLDL) {
375 }
else if (eclDetectorRegion == ECL::DetectorRegion::BWD || eclDetectorRegion == ECL::DetectorRegion::BWDGap) {
376 if (hitStatus == EXT_ECLCROSS) {
378 }
else if (hitStatus == EXT_ECLDL) {
386 return erfc(abs(deltaTheta) / theta_RMS);
391 if (ECL::getDetectorRegion(theta) == ECL::DetectorRegion::BRL) {
392 if (ECL::getDetectorRegion(theta -
m_brlEdgeTheta) != ECL::DetectorRegion::BRL)
return true;
393 else if (ECL::getDetectorRegion(theta +
m_brlEdgeTheta) != ECL::DetectorRegion::BRL)
return true;
400 if (ECL::getDetectorRegion(theta) == ECL::DetectorRegion::FWD || ECL::getDetectorRegion(theta) == ECL::DetectorRegion::FWDGap) {
402 if (pt < matchingThresholdPair.first) {
407 }
else if (ECL::getDetectorRegion(theta) == ECL::DetectorRegion::BRL) {
409 if (theta < matchingThresholdPair.first && pt < matchingThresholdPair.second.first) {
414 }
else if (ECL::getDetectorRegion(theta) == ECL::DetectorRegion::BWD
415 || ECL::getDetectorRegion(theta) == ECL::DetectorRegion::BWDGap) {
417 if (pt < matchingThresholdPair.first) {
static const ChargedStable pion
charged pion particle
@ c_WriteOut
Object/array should be saved by output modules.
@ c_Event
Different object in each event, all objects/arrays are invalidated after event() function has been ca...
Class to store calibrated ECLDigits: ECLCalDigits.
Class to store ECL Showers.
int getHypothesisId() const
Get Hypothesis Id.
double getPhi() const
Get Phi.
int getDetectorRegion() const
Return detector region: 0: below acceptance, 1: FWD, 2: BRL, 3: BWD, 11: FWDGAP, 13: BWDGAP.
double getTheta() const
Get Theta.
TF1 f_thetaRMSBRLNEAR
function to describe theta RMS for BRL NEAR
TF1 f_thetaRMSBRLDL
function to describe theta RMS for BRL DL
TF1 f_thetaRMSBRLCROSS
function to describe theta RMS for BRL CROSS
int m_minimalCDCHits
minimal required number of CDC hits before track-cluster match is initiated
StoreArray< ECLShower > m_eclShowers
Required input array of ECLShowers.
double m_matchingPTThreshold
pt limit between angular-distance based and crystal-entering based matching algorithm
std::vector< std::pair< double, double > > m_matchingThresholdValuesFWD
Matching threshold values for FWD.
std::vector< std::pair< double, std::pair< double, double > > > m_matchingThresholdValuesBRL
Matching threshold values for BRL.
virtual void initialize() override
Use this to initialize resources or memory your module needs.
TF1 f_phiRMSFWDCROSS
function to describe phi RMS for FWD CROSS
DBObjPtr< ECLTrackClusterMatchingParameterizations > m_matchingParameterizations
Parameterizations of RMS.
virtual void event() override
Called once for each event.
StoreArray< TrackFitResult > m_trackFitResults
Required input array of TrackFitResults.
bool m_angularDistanceMatching
members of ECLTrackClusterMatching Module
TF1 f_phiRMSBWDCROSS
function to describe phi RMS for BWD CROSS
virtual void terminate() override
Clean up anything created in initialize().
TF1 f_thetaRMSFWDDL
function to describe theta RMS for FWD DL
double m_matchingConsistency
minimal quality of ExtHit-ECLCluster pair for positive track-cluster match
virtual ~ECLTrackClusterMatchingModule()
Use to clean up anything you created in the constructor.
bool m_useOptimizedMatchingConsistency
if true, a theta dependent matching criterion will be used
TF1 f_phiRMSBWDNEAR
function to describe phi RMS for BWD NEAR
bool isECLHit(const ExtHit &extHit) const
Check if extrapolated hit is inside ECL and matches one of the desired categories.
void optimizedPTMatchingConsistency(double theta, double pt)
choose criterion depending on track's pt
TF1 f_phiRMSBRLDL
function to describe phi RMS for BRL DL
StoreArray< Track > m_tracks
Required input array of Tracks.
bool isECLEnterHit(const ExtHit &extHit) const
Check if status of extrapolated hit is entering of ECL.
TF1 f_phiRMSFWDNEAR
function to describe phi RMS for FWD NEAR
std::vector< std::pair< double, double > > m_matchingThresholdValuesBWD
Matching threshold values for BWD.
TF1 f_thetaRMSBWDDL
function to describe theta RMS for BWD DL
StoreArray< ECLCluster > m_eclClusters
Required input array of ECLClusters.
TF1 f_phiRMSBRLNEAR
function to describe phi RMS for BRL NEAR
TF1 f_thetaRMSBWDCROSS
function to describe theta RMS for BWD CROSS
ECLTrackClusterMatchingModule()
Constructor, for setting module description and parameters.
TF1 f_thetaRMSFWDNEAR
function to describe theta RMS for FWD NEAR
DBObjPtr< ECLTrackClusterMatchingThresholds > m_matchingThresholds
Optimized matching thresholds.
StoreArray< ExtHit > m_extHits
Required input array of ExtHits.
TF1 f_phiRMSFWDDL
function to describe phi RMS for FWD DL
double m_brlEdgeTheta
distance of polar angle from gaps where crystal-entering based matching is applied (in rad)
bool m_skipZeroChargeTracks
if true, tracks whose charge has been set to zero are excluded from track-cluster matching
bool trackTowardsGap(double theta) const
return if track points towards gap or adjacent part of barrel
TF1 f_thetaRMSFWDCROSS
function to describe theta RMS for FWD CROSS
TF1 f_phiRMSBWDDL
function to describe phi RMS for BWD DL
double thetaConsistency(double deltaTheta, double pt, int eclDetectorRegion, int hitStatus) const
Calculate theta consistency based on difference in polar angle.
TF1 f_phiRMSBRLCROSS
function to describe phi RMS for BRL CROSS
TF1 f_thetaRMSBWDNEAR
function to describe theta RMS for BWD NEAR
StoreArray< ECLCalDigit > m_eclCalDigits
Required input array of ECLCalDigits.
double phiConsistency(double deltaPhi, double pt, int eclDetectorRegion, int hitStatus) const
Calculate phi consistency based on difference in azimuthal angle.
double showerQuality(double deltaPhi, double deltaTheta, double pt, int eclDetectorRegion, int hitStatus) const
Calculate matching quality based on phi and theta consistencies.
Store one Ext hit as a ROOT object.
ExtHitStatus getStatus() const
Get state of extrapolation at this hit.
int getCopyID() const
Get detector-element ID of sensitive element within detector.
Const::EDetector getDetectorID() const
Get detector ID of this extrapolation hit.
unsigned short getNHits() const
Get the total Number of CDC hits in the fit.
void setDescription(const std::string &description)
Sets the description of the module.
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
int getArrayIndex() const
Returns this object's array index (in StoreArray), or -1 if not found.
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
FROM * getRelatedFrom(const std::string &name="", const std::string &namedRelation="") const
Get the object from which this object has a relation.
Values of the result of a track fit with a given particle hypothesis.
short getChargeSign() const
Return track charge (1 or -1).
double getTransverseMomentum() const
Getter for the absolute value of the transverse momentum at the perigee.
ROOT::Math::XYZVector getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
Class that bundles various TrackFitResults.
void addParam(const std::string &name, T ¶mVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
ExtHitStatus
Define state of extrapolation for each recorded hit.
Abstract base class for different kinds of events.