Belle II Software development
TrackExtrapolateG4e Class Reference

geant4e-based track extrapolation. More...

#include <TrackExtrapolateG4e.h>

Public Member Functions

 ~TrackExtrapolateG4e ()
 destructor
 
void initialize (double minPt, double minKE, std::vector< Const::ChargedStable > &hypotheses)
 Initialize for track extrapolation by the EXT module.
 
void initialize (double meanDt, double maxDt, double maxSeparation, double maxKLMTrackClusterDistance, double maxECLTrackClusterDistance, double minPt, double minKE, bool addHitsToRecoTrack, std::vector< Const::ChargedStable > &hypotheses)
 Initialize for track extrapolation by the MUID module.
 
void beginRun (bool flag)
 Perform beginning-of-run actions.
 
void event (bool flag)
 Performs track extrapolation for all tracks in one event.
 
void endRun (bool flag)
 Perform end-of-run actions.
 
void terminate (bool flag)
 Terminates this singleton.
 
void extrapolate (int pdgCode, double tof, const G4ThreeVector &position, const G4ThreeVector &momentum, const G4ErrorSymMatrix &covariance)
 Performs track extrapolation for a single track (specified in genfit2 units).
 

Static Public Member Functions

static TrackExtrapolateG4egetInstance ()
 Get the singleton's address.
 

Private Member Functions

 TrackExtrapolateG4e ()
 constructor is hidden; user calls TrackExtrapolateG4e::getInstance() instead
 
 TrackExtrapolateG4e (TrackExtrapolateG4e &)
 copy constructor is hidden; user calls TrackExtrapolateG4e::getInstance() instead
 
void swim (ExtState &, G4ErrorFreeTrajState &, const std::vector< std::pair< ECLCluster *, G4ThreeVector > > *, const std::vector< std::pair< KLMCluster *, G4ThreeVector > > *, std::vector< std::map< const Track *, double > > *)
 Swim a single track (MUID) until it stops or leaves the target cylinder.
 
void swim (ExtState &, G4ErrorFreeTrajState &)
 Swim a single track (EXT) until it stops or leaves the target cylinder.
 
void registerVolumes ()
 Register the list of geant4 physical volumes whose entry/exit points will be saved during extrapolation.
 
void getVolumeID (const G4TouchableHandle &, Const::EDetector &, int &)
 Get the physical volume information for a geant4 physical volume.
 
void fromG4eToPhasespace (const G4ErrorFreeTrajState &, G4ErrorSymMatrix &)
 Convert the geant4e 5x5 covariance to phasespace 6x6 covariance.
 
void fromPhasespaceToG4e (const G4ThreeVector &, const G4ErrorSymMatrix &, G4ErrorTrajErr &)
 Convert the phasespace covariance to geant4e covariance.
 
void fromPhasespaceToG4e (const TVector3 &, const TMatrixDSym &, G4ErrorTrajErr &)
 Convert the phasespace covariance to geant4e covariance.
 
ExtState getStartPoint (const Track &, int, G4ErrorFreeTrajState &)
 Get the start point for a new reconstructed track with specific PDG hypothesis.
 
void createExtHit (ExtHitStatus, const ExtState &, const G4ErrorFreeTrajState &, const G4StepPoint *, const G4TouchableHandle &)
 Create another EXT extrapolation hit for a track candidate.
 
void createECLHit (const ExtState &, const G4ErrorFreeTrajState &, const G4StepPoint *, const G4StepPoint *, const G4TouchableHandle &, const std::pair< ECLCluster *, G4ThreeVector > &, double, double)
 Create another EXT ECL-crystal-crossing hit for a track candidate.
 
bool createMuidHit (ExtState &, G4ErrorFreeTrajState &, KLMMuidLikelihood *, std::vector< std::map< const Track *, double > > *)
 Create another MUID extrapolation hit for a track candidate.
 
bool findBarrelIntersection (ExtState &, const G4ThreeVector &, Intersection &)
 Find the intersection point of the track with the crossed BKLM plane.
 
bool findEndcapIntersection (ExtState &, const G4ThreeVector &, Intersection &)
 Find the intersection point of the track with the crossed EKLM plane.
 
bool findMatchingBarrelHit (Intersection &, const Track *)
 Find the matching BKLM 2D hit nearest the intersection point of the track with the crossed BKLM plane.
 
bool findMatchingEndcapHit (Intersection &, const Track *)
 Find the matching EKLM 2D hit nearest the intersection point of the track with the crossed EKLM plane.
 
void adjustIntersection (Intersection &, const double *, const G4ThreeVector &, const G4ThreeVector &)
 Nudge the track using the matching hit.
 
void finishTrack (const ExtState &, KLMMuidLikelihood *, bool)
 Complete muon identification after end of track extrapolation.
 

Private Attributes

bool m_ExtInitialized
 Flag to indicate that EXT initialize() has been called.
 
bool m_MuidInitialized
 Flag to indicate that MUID initialize() has been called.
 
double m_MeanDt
 Mean hit - trigger time (ns)
 
double m_MaxDt
 Coincidence window half-width for in-time KLM hits (ns)
 
double m_MagneticField
 Magnetic field z component (gauss) at origin.
 
double m_MaxDistSqInVariances
 user-defined maximum squared-distance (in number of variances) for matching hit to extrapolation
 
double m_MaxKLMTrackClusterDistance
 user-defined maximum distance (mm) between KLMCluster and associated track (for KLID)
 
double m_MaxECLTrackClusterDistance
 user-defined maximum distance (mm) between ECLCluster and associated track (for EID)
 
double m_MinPt
 Minimum transverse momentum in MeV/c for extrapolation to be started.
 
double m_MinKE
 Minimum kinetic energy in MeV for extrapolation to continue.
 
Simulation::ExtManagerm_ExtMgr
 Pointer to the ExtManager singleton.
 
const std::vector< Const::ChargedStable > * m_HypothesesExt
 ChargedStable hypotheses for EXT.
 
const std::vector< Const::ChargedStable > * m_HypothesesMuid
 ChargedStable hypotheses for MUID.
 
std::vector< Const::ChargedStable > * m_DefaultHypotheses
 Default ChargedStable hypotheses (needed as call argument but not used)
 
std::map< G4VPhysicalVolume *, enum VolTypes > * m_EnterExit
 Pointers to geant4 physical volumes whose entry/exit points will be saved.
 
std::vector< G4VPhysicalVolume * > * m_BKLMVolumes
 Pointers to BKLM geant4 sensitive (physical) volumes.
 
Simulation::ExtCylSurfaceTargetm_TargetExt
 virtual "target" cylinder for EXT (boundary beyond which extrapolation ends)
 
Simulation::ExtCylSurfaceTargetm_TargetMuid
 virtual "target" cylinder for MUID (boundary beyond which extrapolation ends)
 
DBObjPtr< COILGeometryParm_COILGeometryPar
 Conditions-database object for COIL geometry.
 
DBObjPtr< BeamPipeGeom_BeamPipeGeo
 Conditions-database object for beam pipe geometry.
 
double m_MinRadiusSq
 Minimum squared radius (cm) outside of which extrapolation will continue.
 
double m_OffsetZ
 offset (cm) along z axis of KLM midpoint from IP
 
int m_BarrelNSector
 Number of barrel sectors.
 
double m_BarrelMaxR
 maximum radius (cm) of the barrel
 
double m_BarrelMinR
 minimum radius (cm) of the barrel
 
double m_BarrelHalfLength
 half-length (cm) of the barrel
 
int m_OutermostActiveBarrelLayer
 outermost barrel layer that is active for muon identification (user-defined)
 
double m_BarrelPhiStripVariance [BKLMElementNumbers::getMaximalLayerNumber()+1]
 BKLM RPC phi-measuring strip position variance (cm^2) by layer.
 
double m_BarrelZStripVariance [BKLMElementNumbers::getMaximalLayerNumber()+1]
 BKLM RPC z-measuring strip position variance (cm^2) by layer.
 
double m_BarrelScintVariance
 BKLM scintillator strip position variance (cm^2)
 
double m_BarrelModuleMiddleRadius [2][BKLMElementNumbers::getMaximalSectorNumber()+1][BKLMElementNumbers::getMaximalLayerNumber()+1]
 hit-plane radius (cm) at closest distance to IP of each barrel end | sector | layer
 
G4ThreeVector m_BarrelSectorPerp [BKLMElementNumbers::getMaximalSectorNumber()+1]
 normal unit vector of each barrel sector
 
G4ThreeVector m_BarrelSectorPhi [BKLMElementNumbers::getMaximalSectorNumber()+1]
 azimuthal unit vector of each barrel sector
 
double m_EndcapMaxR
 maximum radius (cm) of the endcaps
 
double m_EndcapMinR
 minimum radius (cm) of the endcaps
 
double m_EndcapMiddleZ
 midpoint along z (cm) of the forward endcap from the KLM midpoint
 
double m_EndcapHalfLength
 half-length (cm) of either endcap
 
int m_OutermostActiveForwardEndcapLayer
 outermost forward-endcap layer that is active for muon identification (user-defined)
 
int m_OutermostActiveBackwardEndcapLayer
 outermost backward-endcap layer that is active for muon identification (user-defined)
 
double m_EndcapScintVariance
 EKLM scintillator strip position variance (cm^2)
 
double m_EndcapModuleMiddleZ [BKLMElementNumbers::getMaximalLayerNumber()+1]
 hit-plane z (cm) of each IP layer relative to KLM midpoint
 
bool m_addHitsToRecoTrack = false
 Parameter to add the found hits also to the reco tracks or not. Is turned off by default.
 
std::map< int, MuidBuilder * > m_MuidBuilderMap
 PDF for the charged final state particle hypotheses.
 
const EKLMElementNumbersm_eklmElementNumbers
 EKLM element numbers.
 
const KLMElementNumbersm_klmElementNumbers
 KLM element numbers.
 
const EKLM::TransformDataGlobalAlignedm_eklmTransformData
 EKLM transformation data.
 
DBObjPtr< KLMChannelStatusm_klmChannelStatus
 Conditions-database object for KLM channel status.
 
DBObjPtr< KLMStripEfficiencym_klmStripEfficiency
 Conditions-database object for KLM strip efficiency.
 
DBObjPtr< KLMLikelihoodParametersm_klmLikelihoodParameters
 Conditions-database object for KLM likelihood parameters.
 
StoreArray< ECLClusterm_eclClusters
 ECL clusters.
 
StoreArray< ExtHitm_extHits
 Ext hits.
 
StoreArray< KLMHit2dm_klmHit2ds
 KLM 2d hits.
 
StoreArray< KLMClusterm_klmClusters
 KLM clusters.
 
StoreArray< KLMMuidHitm_klmMuidHits
 KLM muid hits.
 
StoreArray< KLMMuidLikelihoodm_klmMuidLikelihoods
 KLM muid likelihoods.
 
StoreArray< RecoTrackm_recoTracks
 Reco tracks.
 
StoreArray< Trackm_tracks
 Tracks.
 
StoreArray< TrackClusterSeparationm_trackClusterSeparations
 Track cluster separation.
 

Static Private Attributes

static TrackExtrapolateG4em_Singleton = nullptr
 Stores pointer to the singleton class.
 

Detailed Description

geant4e-based track extrapolation.

This class extrapolates tracks outward from the outer perimeter of the CDC using geant4e.

This class requires a valid geometry in memory (gGeoManager). Therefore, a geometry building module should have been executed before this module is called.

This class has the same functions as a module - and these are called from the ExtModule's so-named functions - but also has an entry that can be called to extrapolate a single user-defined track.

Definition at line 168 of file TrackExtrapolateG4e.h.

Constructor & Destructor Documentation

◆ ~TrackExtrapolateG4e()

destructor

Definition at line 133 of file TrackExtrapolateG4e.cc.

134{
135}

◆ TrackExtrapolateG4e()

TrackExtrapolateG4e ( )
private

constructor is hidden; user calls TrackExtrapolateG4e::getInstance() instead

Definition at line 79 of file TrackExtrapolateG4e.cc.

79 :
80 m_ExtInitialized(false), // initialized later
81 m_MuidInitialized(false), // initialized later
82 m_MeanDt(0.0), // initialized later
83 m_MaxDt(0.0), // initialized later
84 m_MagneticField(0.0), // initialized later
85 m_MaxDistSqInVariances(0.0), // initialized later
86 m_MaxKLMTrackClusterDistance(0.0), // initialized later
87 m_MaxECLTrackClusterDistance(0.0), // initialized later
88 m_MinPt(0.0), // initialized later
89 m_MinKE(0.0), // initialized later
90 m_ExtMgr(nullptr), // initialized later
91 m_HypothesesExt(nullptr), // initialized later
92 m_HypothesesMuid(nullptr), // initialized later
93 m_DefaultHypotheses(nullptr), // initialized later
94 m_EnterExit(nullptr), // initialized later
95 m_BKLMVolumes(nullptr), // initialized later
96 m_TargetExt(nullptr), // initialized later
97 m_TargetMuid(nullptr), // initialized later
98 m_MinRadiusSq(0.0), // initialized later
99 m_OffsetZ(0.0), // initialized later
100 m_BarrelNSector(0), // initialized later
101 m_BarrelMaxR(0.0), // initialized later
102 m_BarrelMinR(0.0), // initialized later
103 m_BarrelHalfLength(0.0), // initialized later
104 m_OutermostActiveBarrelLayer(0), // initialized later
105 m_BarrelScintVariance(0.0), // initialized later
106 m_EndcapMaxR(0.0), // initialized later
107 m_EndcapMinR(0.0), // initialized later
108 m_EndcapMiddleZ(0.0), // initialized later
109 m_EndcapHalfLength(0.0), // initialized later
110 m_OutermostActiveForwardEndcapLayer(0), // initialized later
111 m_OutermostActiveBackwardEndcapLayer(0), // initialized later
112 m_EndcapScintVariance(0.0), // initialized later
113 m_eklmTransformData(nullptr) // initialized later
114{
115 for (int j = 0; j < BKLMElementNumbers::getMaximalLayerNumber() + 1; ++j) {
117 m_BarrelZStripVariance[j] = 0.0;
119 m_EndcapModuleMiddleZ[j] = 0.0;
120 }
121 for (int s = 0; s < BKLMElementNumbers::getMaximalSectorNumber() + 1; ++s) {
122 for (int j = 0; j < BKLMElementNumbers::getMaximalLayerNumber() + 1; ++j) {
123 m_BarrelModuleMiddleRadius[0][s][j] = 0.0;
124 m_BarrelModuleMiddleRadius[1][s][j] = 0.0;
125 }
126 m_BarrelSectorPerp[s] = G4ThreeVector(0.0, 0.0, 0.0);
127 m_BarrelSectorPhi[s] = G4ThreeVector(0.0, 0.0, 0.0);
128 }
131}
static constexpr int getMaximalLayerNumber()
Get maximal layer number (1-based).
static constexpr int getMaximalSectorNumber()
Get maximal sector number (1-based).
static const EKLMElementNumbers & Instance()
Instantiation.
static const KLMElementNumbers & Instance()
Instantiation.
double m_BarrelScintVariance
BKLM scintillator strip position variance (cm^2)
double m_EndcapHalfLength
half-length (cm) of either endcap
double m_MaxKLMTrackClusterDistance
user-defined maximum distance (mm) between KLMCluster and associated track (for KLID)
Simulation::ExtCylSurfaceTarget * m_TargetMuid
virtual "target" cylinder for MUID (boundary beyond which extrapolation ends)
int m_OutermostActiveForwardEndcapLayer
outermost forward-endcap layer that is active for muon identification (user-defined)
double m_BarrelMaxR
maximum radius (cm) of the barrel
G4ThreeVector m_BarrelSectorPhi[BKLMElementNumbers::getMaximalSectorNumber()+1]
azimuthal unit vector of each barrel sector
int m_OutermostActiveBackwardEndcapLayer
outermost backward-endcap layer that is active for muon identification (user-defined)
std::vector< Const::ChargedStable > * m_DefaultHypotheses
Default ChargedStable hypotheses (needed as call argument but not used)
const std::vector< Const::ChargedStable > * m_HypothesesExt
ChargedStable hypotheses for EXT.
double m_MagneticField
Magnetic field z component (gauss) at origin.
double m_BarrelHalfLength
half-length (cm) of the barrel
double m_MaxDt
Coincidence window half-width for in-time KLM hits (ns)
double m_BarrelMinR
minimum radius (cm) of the barrel
double m_BarrelModuleMiddleRadius[2][BKLMElementNumbers::getMaximalSectorNumber()+1][BKLMElementNumbers::getMaximalLayerNumber()+1]
hit-plane radius (cm) at closest distance to IP of each barrel end | sector | layer
double m_BarrelPhiStripVariance[BKLMElementNumbers::getMaximalLayerNumber()+1]
BKLM RPC phi-measuring strip position variance (cm^2) by layer.
double m_MinRadiusSq
Minimum squared radius (cm) outside of which extrapolation will continue.
double m_EndcapMiddleZ
midpoint along z (cm) of the forward endcap from the KLM midpoint
const EKLMElementNumbers * m_eklmElementNumbers
EKLM element numbers.
const std::vector< Const::ChargedStable > * m_HypothesesMuid
ChargedStable hypotheses for MUID.
double m_OffsetZ
offset (cm) along z axis of KLM midpoint from IP
const KLMElementNumbers * m_klmElementNumbers
KLM element numbers.
std::vector< G4VPhysicalVolume * > * m_BKLMVolumes
Pointers to BKLM geant4 sensitive (physical) volumes.
int m_OutermostActiveBarrelLayer
outermost barrel layer that is active for muon identification (user-defined)
G4ThreeVector m_BarrelSectorPerp[BKLMElementNumbers::getMaximalSectorNumber()+1]
normal unit vector of each barrel sector
double m_EndcapModuleMiddleZ[BKLMElementNumbers::getMaximalLayerNumber()+1]
hit-plane z (cm) of each IP layer relative to KLM midpoint
double m_MaxECLTrackClusterDistance
user-defined maximum distance (mm) between ECLCluster and associated track (for EID)
double m_EndcapScintVariance
EKLM scintillator strip position variance (cm^2)
double m_MinPt
Minimum transverse momentum in MeV/c for extrapolation to be started.
int m_BarrelNSector
Number of barrel sectors.
std::map< G4VPhysicalVolume *, enum VolTypes > * m_EnterExit
Pointers to geant4 physical volumes whose entry/exit points will be saved.
bool m_MuidInitialized
Flag to indicate that MUID initialize() has been called.
double m_EndcapMinR
minimum radius (cm) of the endcaps
double m_MeanDt
Mean hit - trigger time (ns)
double m_EndcapMaxR
maximum radius (cm) of the endcaps
Simulation::ExtCylSurfaceTarget * m_TargetExt
virtual "target" cylinder for EXT (boundary beyond which extrapolation ends)
bool m_ExtInitialized
Flag to indicate that EXT initialize() has been called.
double m_MinKE
Minimum kinetic energy in MeV for extrapolation to continue.
Simulation::ExtManager * m_ExtMgr
Pointer to the ExtManager singleton.
double m_MaxDistSqInVariances
user-defined maximum squared-distance (in number of variances) for matching hit to extrapolation
double m_BarrelZStripVariance[BKLMElementNumbers::getMaximalLayerNumber()+1]
BKLM RPC z-measuring strip position variance (cm^2) by layer.
const EKLM::TransformDataGlobalAligned * m_eklmTransformData
EKLM transformation data.

Member Function Documentation

◆ adjustIntersection()

void adjustIntersection ( Intersection ,
const double *  ,
const G4ThreeVector &  ,
const G4ThreeVector &   
)
private

Nudge the track using the matching hit.

Definition at line 1639 of file TrackExtrapolateG4e.cc.

1641{
1642 // Use the gain matrix formalism to get the corrected track parameters.
1643 // R. Fruhwirth, Application of Kalman Filtering, NIM A262 (1987) 444-450
1644 // Equations (7)
1645 // x_k^{k-1} = extPar[] 6 elements before filtering
1646 // C_k^{k-1} = extCov[] 6x6 symmetric before filtering
1647 // r_k^{k-1} = residual[] 2 elements before filtering
1648 // h_k = 2x6 projects cartesian coordinates to measurement-plane coordinates
1649 // H_k = @h_k/@x = jacobian[] 2x6 Jacobian of projection to measurement plane
1650 // R_k^{k-1} = correction[] 2x2 before Invert()
1651 // G_k = R^(-1) = correction[] 2x2 after Invert()
1652 // K_k = gain[] 6x2 Kalman gain matrix
1653 // x_k = extPar[] 6 elements after filtering
1654 // C_k = extCov[] 6x6 symmetric after filtering
1655 // r_k = residual[] 2 elements after filtering
1656 // Use the relation K*H*C = (C*H^T*R^-1)*H*C = C*(H^T*R^-1*H)*C^T
1657
1658 // In most cases, extPos0 is the same as intersection.position. They differ only when
1659 // the nearest BKLM hit is in the sector adjacent to that of intersection.position.
1660 G4ThreeVector extPos(extPos0);
1661 G4ThreeVector extMom(intersection.momentum);
1662 G4ThreeVector extDir(extMom.unit());
1663 G4ThreeVector diffPos(hitPos - extPos);
1664 G4ErrorSymMatrix extCov(intersection.covariance);
1665 // Track parameters (x,y,z,px,py,pz) before correction
1666 G4ErrorMatrix extPar(6, 1); // initialized to all zeroes
1667 extPar[0][0] = extPos.x();
1668 extPar[1][0] = extPos.y();
1669 extPar[2][0] = extPos.z();
1670 extPar[3][0] = extMom.x();
1671 extPar[4][0] = extMom.y();
1672 extPar[5][0] = extMom.z();
1673 G4ThreeVector nA; // unit vector normal to the readout plane
1674 G4ThreeVector nB; // unit vector along phi- or x-readout direction (for barrel or endcap)
1675 G4ThreeVector nC; // unit vector along z- or y-readout direction (for barrel or endcap)
1676 if (intersection.inBarrel) {
1677 nA = m_BarrelSectorPerp[intersection.sector];
1678 nB = m_BarrelSectorPhi[intersection.sector];
1679 nC = G4ThreeVector(0.0, 0.0, 1.0);
1680 } else {
1681 double out = (intersection.isForward ? 1.0 : -1.0);
1682 nA = G4ThreeVector(0.0, 0.0, out);
1683 nB = G4ThreeVector(out, 0.0, 0.0);
1684 nC = G4ThreeVector(0.0, out, 0.0);
1685 }
1686 // Don't adjust the extrapolation if the track is nearly tangent to the readout plane.
1687 double extDirA = extDir * nA;
1688 if (std::fabs(extDirA) < 1.0E-6)
1689 return;
1690 double extDirBA = extDir * nB / extDirA;
1691 double extDirCA = extDir * nC / extDirA;
1692 // Move the extrapolated coordinate (at most a tiny amount!) to the plane of the hit.
1693 // If the moved point is outside the KLM, don't do Kalman filtering.
1694 G4ThreeVector move = extDir * ((diffPos * nA) / extDirA);
1695 extPos += move;
1696 diffPos -= move;
1697 intersection.positionAtHitPlane = extPos;
1698 // Projection jacobian onto the nB-nC measurement plane
1699 G4ErrorMatrix jacobian(2, 6); // initialized to all zeroes
1700 jacobian[0][0] = nB.x() - nA.x() * extDirBA;
1701 jacobian[0][1] = nB.y() - nA.y() * extDirBA;
1702 jacobian[0][2] = nB.z() - nA.z() * extDirBA;
1703 jacobian[1][0] = nC.x() - nA.x() * extDirCA;
1704 jacobian[1][1] = nC.y() - nA.y() * extDirCA;
1705 jacobian[1][2] = nC.z() - nA.z() * extDirCA;
1706 // Residuals of EXT track and KLM hit on the nB-nC measurement plane
1707 G4ErrorMatrix residual(2, 1); // initialized to all zeroes
1708 residual[0][0] = diffPos.x() * jacobian[0][0] + diffPos.y() * jacobian[0][1] + diffPos.z() * jacobian[0][2];
1709 residual[1][0] = diffPos.x() * jacobian[1][0] + diffPos.y() * jacobian[1][1] + diffPos.z() * jacobian[1][2];
1710 // Measurement errors in the detector plane
1711 G4ErrorSymMatrix hitCov(2, 0); // initialized to all zeroes
1712 hitCov[0][0] = localVariance[0];
1713 hitCov[1][1] = localVariance[1];
1714 // No magnetic field: increase the hit uncertainty
1715 if (m_MagneticField == 0.0) {
1716 hitCov[0][0] *= 10.0;
1717 hitCov[1][1] *= 10.0;
1718 }
1719 // Now get the correction matrix: combined covariance of EXT and KLM hit.
1720 // 1st dimension = nB, 2nd dimension = nC.
1721 G4ErrorSymMatrix correction(extCov.similarity(jacobian) + hitCov);
1722 // Ignore the best hit if it is too far from the extrapolated-track intersection in the hit's plane
1723 if (residual[0][0] * residual[0][0] > correction[0][0] * m_MaxDistSqInVariances)
1724 return;
1725 if (residual[1][0] * residual[1][0] > correction[1][1] * m_MaxDistSqInVariances)
1726 return;
1727 int fail = 0;
1728 correction.invert(fail);
1729 if (fail != 0)
1730 return;
1731 // Matrix inversion succeeded and is reasonable.
1732 // Evaluate chi-squared increment assuming that the Kalman filter
1733 // won't be able to adjust the extrapolated track's position (fall-back).
1734 intersection.chi2 = (correction.similarityT(residual))[0][0];
1735 // Do the Kalman filtering
1736 G4ErrorMatrix gain((extCov * jacobian.T()) * correction);
1737 G4ErrorSymMatrix HRH(correction.similarityT(jacobian));
1738 extCov -= HRH.similarity(extCov);
1739 extPar += gain * residual;
1740 extPos.set(extPar[0][0], extPar[1][0], extPar[2][0]);
1741 extMom.set(extPar[3][0], extPar[4][0], extPar[5][0]);
1742 // Calculate the chi-squared increment using the Kalman-filtered state
1743 correction = hitCov - extCov.similarity(jacobian);
1744 correction.invert(fail);
1745 if (fail != 0)
1746 return;
1747 diffPos = hitPos - extPos;
1748 residual[0][0] = diffPos.x() * jacobian[0][0] + diffPos.y() * jacobian[0][1] + diffPos.z() * jacobian[0][2];
1749 residual[1][0] = diffPos.x() * jacobian[1][0] + diffPos.y() * jacobian[1][1] + diffPos.z() * jacobian[1][2];
1750 intersection.chi2 = (correction.similarityT(residual))[0][0];
1751 // Update the position, momentum and covariance of the point
1752 // Project the corrected extrapolation to the plane of the original
1753 // extrapolation's intersection.position. (Note: intersection.position is the same as
1754 // extPos0 in all cases except when nearest BKLM hit is in adjacent
1755 // sector, for which extPos0 is a projected position to the hit's plane.)
1756 // Also, leave the momentum magnitude unchanged.
1757 intersection.position = extPos + extDir * (((intersection.position - extPos) * nA) / extDirA);
1758 intersection.momentum = intersection.momentum.mag() * extMom.unit();
1759 intersection.covariance = extCov;
1760}

◆ beginRun()

void beginRun ( bool  flag)

Perform beginning-of-run actions.

Parameters
flagTrue if called by Muid module, false if called by Ext module.

Definition at line 333 of file TrackExtrapolateG4e.cc.

334{
335 B2DEBUG(20, (byMuid ? "muid" : "ext"));
336 if (byMuid) {
337 if (!m_klmChannelStatus.isValid())
338 B2FATAL("KLM channel status data are not available.");
339 if (!m_klmStripEfficiency.isValid())
340 B2FATAL("KLM strip efficiency data are not available.");
341 if (!m_klmLikelihoodParameters.isValid())
342 B2FATAL("KLM likelihood parameters are not available.");
343 std::vector<int> muidPdgCodes = MuidElementNumbers::getPDGVector();
344 if (!m_MuidBuilderMap.empty()) {
345 if (m_klmLikelihoodParameters.hasChanged()) { /* Clear m_MuidBuilderMap if KLMLikelihoodParameters payload changed. */
346 for (auto const& muidBuilder : m_MuidBuilderMap)
347 delete muidBuilder.second;
348 m_MuidBuilderMap.clear();
349 } else /* Return if m_MuidBuilderMap is already initialized. */
350 return;
351 }
352 for (int pdg : muidPdgCodes)
353 m_MuidBuilderMap.insert(std::pair<int, MuidBuilder*>(pdg, new MuidBuilder(pdg)));
354 }
355}
Build the Muid likelihoods starting from the hit pattern and the transverse scattering in the KLM.
Definition: MuidBuilder.h:29
static std::vector< int > getPDGVector()
Get a vector with all the hypothesis PDG codes used for Muid.
DBObjPtr< KLMLikelihoodParameters > m_klmLikelihoodParameters
Conditions-database object for KLM likelihood parameters.
std::map< int, MuidBuilder * > m_MuidBuilderMap
PDF for the charged final state particle hypotheses.
DBObjPtr< KLMChannelStatus > m_klmChannelStatus
Conditions-database object for KLM channel status.
DBObjPtr< KLMStripEfficiency > m_klmStripEfficiency
Conditions-database object for KLM strip efficiency.

◆ createExtHit()

void createExtHit ( ExtHitStatus  status,
const ExtState extState,
const G4ErrorFreeTrajState &  g4eState,
const G4StepPoint *  stepPoint,
const G4TouchableHandle &  touch 
)
private

Create another EXT extrapolation hit for a track candidate.

Definition at line 1182 of file TrackExtrapolateG4e.cc.

1185{
1186 Const::EDetector detID(Const::EDetector::invalidDetector);
1187 int copyID(0);
1188 getVolumeID(touch, detID, copyID);
1189 G4ThreeVector pos(stepPoint->GetPosition() / CLHEP::cm);
1190 G4ThreeVector mom(stepPoint->GetMomentum() / CLHEP::GeV);
1191 if (extState.isCosmic)
1192 mom = -mom;
1193 G4ErrorSymMatrix covariance(6, 0);
1194 fromG4eToPhasespace(g4eState, covariance);
1195 ExtHit* extHit = m_extHits.appendNew(extState.pdgCode, detID, copyID, status,
1196 extState.isCosmic, extState.tof,
1197 pos, mom, covariance);
1198 // If called standalone, there will be no associated track
1199 if (extState.track != nullptr)
1200 extState.track->addRelationTo(extHit);
1201}
EDetector
Enum for identifying the detector components (detector and subdetector).
Definition: Const.h:42
Store one Ext hit as a ROOT object.
Definition: ExtHit.h:31
void addRelationTo(const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
Add a relation from this object to another object (with caching).
void getVolumeID(const G4TouchableHandle &, Const::EDetector &, int &)
Get the physical volume information for a geant4 physical volume.
StoreArray< ExtHit > m_extHits
Ext hits.
void fromG4eToPhasespace(const G4ErrorFreeTrajState &, G4ErrorSymMatrix &)
Convert the geant4e 5x5 covariance to phasespace 6x6 covariance.
bool isCosmic
True for back-propagation of a cosmic ray.
int pdgCode
Particle hypothesis that is being extrapolated.
double tof
Time of flight from IP (ns), updated during extrapolation.
const Track * track
Pointer to the reconstructed track.

◆ createMuidHit()

bool createMuidHit ( ExtState extState,
G4ErrorFreeTrajState &  g4eState,
KLMMuidLikelihood klmMuidLikelihood,
std::vector< std::map< const Track *, double > > *  bklmHitUsed 
)
private

Create another MUID extrapolation hit for a track candidate.

Definition at line 1206 of file TrackExtrapolateG4e.cc.

1208{
1209
1210 Intersection intersection;
1211 intersection.hit = -1;
1212 intersection.chi2 = -1.0;
1213 intersection.position = g4eState.GetPosition() / CLHEP::cm;
1214 intersection.momentum = g4eState.GetMomentum() / CLHEP::GeV;
1215 G4ThreeVector prePos = g4eState.GetG4Track()->GetStep()->GetPreStepPoint()->GetPosition() / CLHEP::cm;
1216 G4ThreeVector oldPosition(prePos.x(), prePos.y(), prePos.z());
1217 double r = intersection.position.perp();
1218 double z = std::fabs(intersection.position.z() - m_OffsetZ);
1219
1220 // Is the track in the barrel?
1221 if ((r > m_BarrelMinR) && (r < m_BarrelMaxR) && (z < m_BarrelHalfLength)) {
1222 // Did the track cross the inner midplane of a detector module?
1223 if (findBarrelIntersection(extState, oldPosition, intersection)) {
1224 fromG4eToPhasespace(g4eState, intersection.covariance);
1225 if (findMatchingBarrelHit(intersection, extState.track)) {
1226 (*bklmHitUsed)[intersection.hit].insert(std::pair<const Track*, double>(extState.track, intersection.chi2));
1227 extState.extLayerPattern |= (0x00000001 << intersection.layer);
1228 float layerBarrelEfficiency = 1.;
1229 for (int plane = 0; plane <= BKLMElementNumbers::getMaximalPlaneNumber(); plane++) {
1230 layerBarrelEfficiency *= m_klmStripEfficiency->getBarrelEfficiency(
1232 intersection.sector + 1, intersection.layer + 1, plane, 1);
1233 }
1234 klmMuidLikelihood->setExtBKLMEfficiencyValue(intersection.layer, layerBarrelEfficiency);
1235 if (extState.lastBarrelExtLayer < intersection.layer) {
1236 extState.lastBarrelExtLayer = intersection.layer;
1237 }
1238 extState.hitLayerPattern |= (0x00000001 << intersection.layer);
1239 if (extState.lastBarrelHitLayer < intersection.layer) {
1240 extState.lastBarrelHitLayer = intersection.layer;
1241 }
1242 // If the updated point is outside the barrel, discard it and the Kalman-fitter adjustment
1243 r = intersection.position.perp();
1244 z = std::fabs(intersection.position.z() - m_OffsetZ);
1245 if ((r <= m_BarrelMinR) || (r >= m_BarrelMaxR) || (z >= m_BarrelHalfLength)) {
1246 intersection.chi2 = -1.0;
1247 }
1248 } else {
1249 // Record a no-hit track crossing if this step is strictly within a barrel sensitive volume
1250 std::vector<G4VPhysicalVolume*>::iterator j = find(m_BKLMVolumes->begin(), m_BKLMVolumes->end(),
1251 g4eState.GetG4Track()->GetVolume());
1252 if (j != m_BKLMVolumes->end()) {
1253 bool isDead = true; // by default, the nearest orthogonal strips are dead
1254 int section = intersection.isForward ?
1257 int sector = intersection.sector + 1; // from 0-based to 1-based enumeration
1258 int layer = intersection.layer + 1; // from 0-based to 1-based enumeration
1259 const bklm::Module* m = bklm::GeometryPar::instance()->findModule(section, sector, layer); // uses 1-based enumeration
1260 if (m) {
1261 const CLHEP::Hep3Vector localPosition = m->globalToLocal(intersection.position); // uses and returns position in cm
1262 int zStrip = m->getZStripNumber(localPosition);
1263 int phiStrip = m->getPhiStripNumber(localPosition);
1264 if (zStrip >= 0 && phiStrip >= 0) {
1265 uint16_t channel1, channel2;
1267 section, sector, layer,
1270 section, sector, layer,
1272 enum KLMChannelStatus::ChannelStatus status1, status2;
1273 status1 = m_klmChannelStatus->getChannelStatus(channel1);
1274 status2 = m_klmChannelStatus->getChannelStatus(channel2);
1275 if (status1 == KLMChannelStatus::c_Unknown ||
1276 status2 == KLMChannelStatus::c_Unknown)
1277 B2ERROR("No KLM channel status data."
1278 << LogVar("Section", section)
1279 << LogVar("Sector", sector)
1280 << LogVar("Layer", layer)
1281 << LogVar("Z strip", zStrip)
1282 << LogVar("Phi strip", phiStrip));
1283 isDead = (status1 == KLMChannelStatus::c_Dead ||
1284 status2 == KLMChannelStatus::c_Dead);
1285 }
1286 }
1287 if (!isDead) {
1288 extState.extLayerPattern |= (0x00000001 << intersection.layer); // valid extrapolation-crossing of the layer but no matching hit
1289 float layerBarrelEfficiency = 1.;
1290 for (int plane = 0; plane <= BKLMElementNumbers::getMaximalPlaneNumber(); plane++) {
1291 layerBarrelEfficiency *= m_klmStripEfficiency->getBarrelEfficiency(
1293 intersection.sector + 1, intersection.layer + 1, plane, 1);
1294 }
1295 klmMuidLikelihood->setExtBKLMEfficiencyValue(intersection.layer, layerBarrelEfficiency);
1296 } else {
1297 klmMuidLikelihood->setExtBKLMEfficiencyValue(intersection.layer, 0);
1298 }
1299 if (extState.lastBarrelExtLayer < intersection.layer) {
1300 extState.lastBarrelExtLayer = intersection.layer;
1301 }
1302 }
1303 }
1304 }
1305 }
1306
1307 // Is the track in the endcap?
1308 if ((r > m_EndcapMinR) && (std::fabs(z - m_EndcapMiddleZ) < m_EndcapHalfLength)) {
1309 // Did the track cross the inner midplane of a detector module?
1310 if (findEndcapIntersection(extState, oldPosition, intersection)) {
1311 fromG4eToPhasespace(g4eState, intersection.covariance);
1312 if (findMatchingEndcapHit(intersection, extState.track)) {
1313 extState.extLayerPattern |= (0x00008000 << intersection.layer);
1314 float layerEndcapEfficiency = 1.;
1315 for (int plane = 1; plane <= EKLMElementNumbers::getMaximalPlaneNumber(); plane++) {
1316 layerEndcapEfficiency *= m_klmStripEfficiency->getEndcapEfficiency(
1318 intersection.sector + 1, intersection.layer + 1, plane, 1);
1319 }
1320 klmMuidLikelihood->setExtEKLMEfficiencyValue(intersection.layer, layerEndcapEfficiency);
1321 if (extState.lastEndcapExtLayer < intersection.layer) {
1322 extState.lastEndcapExtLayer = intersection.layer;
1323 }
1324 extState.hitLayerPattern |= (0x00008000 << intersection.layer);
1325 if (extState.lastEndcapHitLayer < intersection.layer) {
1326 extState.lastEndcapHitLayer = intersection.layer;
1327 }
1328 // If the updated point is outside the endcap, discard it and the Kalman-fitter adjustment
1329 r = intersection.position.perp();
1330 z = std::fabs(intersection.position.z() - m_OffsetZ);
1331 if ((r <= m_EndcapMinR) || (r >= m_EndcapMaxR) || (std::fabs(z - m_EndcapMiddleZ) >= m_EndcapHalfLength)) {
1332 intersection.chi2 = -1.0;
1333 }
1334 } else {
1335 bool isDead = true;
1336 int result, strip1, strip2;
1338 intersection.position, &strip1, &strip2);
1339 if (result == 0) {
1340 uint16_t channel1, channel2;
1341 channel1 = m_klmElementNumbers->channelNumberEKLM(strip1);
1342 channel2 = m_klmElementNumbers->channelNumberEKLM(strip2);
1343 enum KLMChannelStatus::ChannelStatus status1, status2;
1344 status1 = m_klmChannelStatus->getChannelStatus(channel1);
1345 status2 = m_klmChannelStatus->getChannelStatus(channel2);
1346 if (status1 == KLMChannelStatus::c_Unknown ||
1347 status2 == KLMChannelStatus::c_Unknown)
1348 B2ERROR("Incomplete KLM channel status data.");
1349 isDead = (status1 == KLMChannelStatus::c_Dead ||
1350 status2 == KLMChannelStatus::c_Dead);
1351 }
1352 if (!isDead) {
1353 extState.extLayerPattern |= (0x00008000 << intersection.layer); // valid extrapolation-crossing of the layer but no matching hit
1354 float layerEndcapEfficiency = 1.;
1355 for (int plane = 1; plane <= EKLMElementNumbers::getMaximalPlaneNumber(); plane++) {
1356 layerEndcapEfficiency *= m_klmStripEfficiency->getEndcapEfficiency(
1358 intersection.sector + 1, intersection.layer + 1, plane, 1);
1359 }
1360 klmMuidLikelihood->setExtEKLMEfficiencyValue(intersection.layer, layerEndcapEfficiency);
1361 } else {
1362 klmMuidLikelihood->setExtEKLMEfficiencyValue(intersection.layer, 0);
1363 }
1364 if (extState.lastEndcapExtLayer < intersection.layer) {
1365 extState.lastEndcapExtLayer = intersection.layer;
1366 }
1367 }
1368 }
1369 }
1370
1371 // Create a new MuidHit and RelationEntry between it and the track.
1372 // Adjust geant4e's position, momentum and covariance based on matching hit and tell caller to update the geant4e state.
1373 if (intersection.chi2 >= 0.0) {
1374 ROOT::Math::XYZVector tpos(intersection.position.x(),
1375 intersection.position.y(),
1376 intersection.position.z());
1377 ROOT::Math::XYZVector tposAtHitPlane(intersection.positionAtHitPlane.x(),
1378 intersection.positionAtHitPlane.y(),
1379 intersection.positionAtHitPlane.z());
1380 KLMMuidHit* klmMuidHit = m_klmMuidHits.appendNew(extState.pdgCode,
1381 intersection.inBarrel, intersection.isForward,
1382 intersection.sector, intersection.layer,
1383 tpos, tposAtHitPlane,
1384 extState.tof, intersection.time,
1385 intersection.chi2);
1386 if (extState.track != nullptr) { extState.track->addRelationTo(klmMuidHit); }
1387 G4Point3D newPos(intersection.position.x() * CLHEP::cm,
1388 intersection.position.y() * CLHEP::cm,
1389 intersection.position.z() * CLHEP::cm);
1390 g4eState.SetPosition(newPos);
1391 G4Vector3D newMom(intersection.momentum.x() * CLHEP::GeV,
1392 intersection.momentum.y() * CLHEP::GeV,
1393 intersection.momentum.z() * CLHEP::GeV);
1394 g4eState.SetMomentum(newMom);
1395 G4ErrorTrajErr covG4e;
1396 fromPhasespaceToG4e(intersection.momentum, intersection.covariance, covG4e);
1397 g4eState.SetError(covG4e);
1398 extState.chi2 += intersection.chi2;
1399 extState.nPoint += 2; // two (orthogonal) independent hits per detector layer
1400 return true;
1401 }
1402
1403 // Tell caller that the geant4e state was not modified.
1404 return false;
1405
1406}
static constexpr int getMaximalPlaneNumber()
Get maximal plane number (0-based).
static constexpr int getMaximalPlaneNumber()
Get maximal plane number.
int getStripsByIntersection(const HepGeom::Point3D< double > &intersection, int *strip1, int *strip2) const
Find strips by intersection.
ChannelStatus
Channel status.
@ c_Dead
Dead channel (no signal).
@ c_Unknown
Unknown status (no data).
KLMChannelNumber channelNumberBKLM(int section, int sector, int layer, int plane, int strip) const
Get channel number for BKLM.
KLMChannelNumber channelNumberEKLM(int section, int sector, int layer, int plane, int strip) const
Get channel number for EKLM.
Store one muon-identification hit in the KLM as a ROOT object.
Definition: KLMMuidHit.h:24
void setExtEKLMEfficiencyValue(int layer, float efficiency)
Set the efficiency of a given EKLM layer.
void setExtBKLMEfficiencyValue(int layer, float efficiency)
Set the efficiency of a given BKLM layer.
bool findEndcapIntersection(ExtState &, const G4ThreeVector &, Intersection &)
Find the intersection point of the track with the crossed EKLM plane.
StoreArray< KLMMuidHit > m_klmMuidHits
KLM muid hits.
bool findMatchingEndcapHit(Intersection &, const Track *)
Find the matching EKLM 2D hit nearest the intersection point of the track with the crossed EKLM plane...
void fromPhasespaceToG4e(const G4ThreeVector &, const G4ErrorSymMatrix &, G4ErrorTrajErr &)
Convert the phasespace covariance to geant4e covariance.
bool findMatchingBarrelHit(Intersection &, const Track *)
Find the matching BKLM 2D hit nearest the intersection point of the track with the crossed BKLM plane...
bool findBarrelIntersection(ExtState &, const G4ThreeVector &, Intersection &)
Find the intersection point of the track with the crossed BKLM plane.
const Module * findModule(int section, int sector, int layer) const
Get the pointer to the definition of a module.
Definition: GeometryPar.cc:721
static GeometryPar * instance(void)
Static method to get a reference to the singleton GeometryPar instance.
Definition: GeometryPar.cc:27
Define the geometry of a BKLM module Each sector [octant] contains Modules.
Definition: Module.h:76
Class to store variables with their name which were sent to the logging service.
int lastEndcapExtLayer
MUID: outermost endcap layer crossed by the extrapolated track.
double chi2
MUID: accumulated chi-squared of all in-plane transverse deviations between extrapolation and matchin...
int nPoint
MUID: accumulated number of points with matching 2D hits.
int lastEndcapHitLayer
MUID: outermost endcap layer with a matching hit.
int extLayerPattern
MUID: accumulated bit pattern of layers crossed by the extrapolated track.
int lastBarrelHitLayer
MUID: outermost barrel layer with a matching hit.
int hitLayerPattern
MUID: accumulated bit pattern of layers with matching hits.
int lastBarrelExtLayer
MUID: outermost barrel layer crossed by the extrapolated track.
intersection of muid-extrapolated track with a KLM layer
int sector
sector number (0..7 for barrel, 0..3 for endcap) of this point
G4ThreeVector momentum
extrapolated-track momentum (GeV/c) at this intersection
double chi2
chi-squared value of transverse deviation between extrapolated and measured hit positions
bool inBarrel
flag to indicate if this point is in the barrel (true) or endcap (false)
int hit
index in {B,E}KLMHit2ds of matching hit
G4ThreeVector positionAtHitPlane
extrapolated-track position (cm) projected to the 2D hit's midplane
G4ThreeVector position
extrapolated-track global position (cm) of this intersection
G4ErrorSymMatrix covariance
extrapolated-track phase-space covariance matrix at this intersection
bool isForward
flag to indicate if this point is in the forward (true) or backward (false) end
double time
time (ns) of matching BKLMHit2d
int layer
layer number (0..14 for barrel, 0..13 for endcap) of this point

◆ endRun()

void endRun ( bool  flag)

Perform end-of-run actions.

Parameters
flagTrue if called by Muid module, false if called by Ext module.

Definition at line 413 of file TrackExtrapolateG4e.cc.

414{
415}

◆ event()

void event ( bool  flag)

Performs track extrapolation for all tracks in one event.

Parameters
flagTrue if called by Muid module, false if called by Ext module.

Definition at line 357 of file TrackExtrapolateG4e.cc.

358{
359
360 // Put geant4 in proper state (in case this module is in a separate process)
361 if (G4StateManager::GetStateManager()->GetCurrentState() == G4State_Idle) {
362 G4StateManager::GetStateManager()->SetNewState(G4State_GeomClosed);
363 }
364
365 G4ThreeVector directionAtIP, positionG4e, momentumG4e;
366 G4ErrorTrajErr covG4e(5); // initialized to zeroes
367
368 // Loop over the reconstructed tracks
369 // Do extrapolation for each hypothesis of each reconstructed track.
370 if (byMuid) { // event() called by Muid module
371 G4ErrorPropagatorData::GetErrorPropagatorData()->SetTarget(m_TargetMuid);
372 std::vector<std::pair<ECLCluster*, G4ThreeVector> > eclClusterInfo(m_eclClusters.getEntries());
373 for (int c = 0; c < m_eclClusters.getEntries(); ++c) {
374 eclClusterInfo[c].first = m_eclClusters[c];
375 eclClusterInfo[c].second = G4ThreeVector(m_eclClusters[c]->getClusterPosition().X(),
376 m_eclClusters[c]->getClusterPosition().Y(),
377 m_eclClusters[c]->getClusterPosition().Z()) * CLHEP::cm;
378 }
379 std::vector<std::pair<KLMCluster*, G4ThreeVector> > klmClusterInfo(m_klmClusters.getEntries());
380 for (int c = 0; c < m_klmClusters.getEntries(); ++c) {
381 klmClusterInfo[c].first = m_klmClusters[c];
382 klmClusterInfo[c].second = G4ThreeVector(m_klmClusters[c]->getClusterPosition().X(),
383 m_klmClusters[c]->getClusterPosition().Y(),
384 m_klmClusters[c]->getClusterPosition().Z()) * CLHEP::cm;
385 }
386 // Keep track of (re-)use of BKLMHit2ds
387 std::vector<std::map<const Track*, double> > bklmHitUsed(m_klmHit2ds.getEntries());
388 for (auto& b2track : m_tracks) {
389 for (const auto& hypothesis : *m_HypothesesMuid) {
390 int pdgCode = hypothesis.getPDGCode();
391 if (hypothesis == Const::electron || hypothesis == Const::muon)
392 pdgCode = -pdgCode;
393 G4ErrorFreeTrajState g4eState("g4e_mu+", G4ThreeVector(), G4ThreeVector()); // will be updated
394 ExtState extState = getStartPoint(b2track, pdgCode, g4eState);
395 swim(extState, g4eState, &eclClusterInfo, &klmClusterInfo, &bklmHitUsed);
396 } // Muid hypothesis loop
397 } // Muid track loop
398 } else { // event() called by Ext module
399 G4ErrorPropagatorData::GetErrorPropagatorData()->SetTarget(m_TargetExt);
400 for (auto& b2track : m_tracks) {
401 for (const auto& hypothesis : *m_HypothesesExt) {
402 int pdgCode = hypothesis.getPDGCode();
403 if (hypothesis == Const::electron || hypothesis == Const::muon) pdgCode = -pdgCode;
404 G4ErrorFreeTrajState g4eState("g4e_mu+", G4ThreeVector(), G4ThreeVector()); // will be updated
405 ExtState extState = getStartPoint(b2track, pdgCode, g4eState);
406 swim(extState, g4eState);
407 } // Ext hypothesis loop
408 } // Ext track loop
409 } // byMuid
410
411}
static const ChargedStable muon
muon particle
Definition: Const.h:660
static const ChargedStable electron
electron particle
Definition: Const.h:659
ExtState getStartPoint(const Track &, int, G4ErrorFreeTrajState &)
Get the start point for a new reconstructed track with specific PDG hypothesis.
StoreArray< KLMCluster > m_klmClusters
KLM clusters.
StoreArray< Track > m_tracks
Tracks.
StoreArray< ECLCluster > m_eclClusters
ECL clusters.
StoreArray< KLMHit2d > m_klmHit2ds
KLM 2d hits.
void swim(ExtState &, G4ErrorFreeTrajState &, const std::vector< std::pair< ECLCluster *, G4ThreeVector > > *, const std::vector< std::pair< KLMCluster *, G4ThreeVector > > *, std::vector< std::map< const Track *, double > > *)
Swim a single track (MUID) until it stops or leaves the target cylinder.
Data structure to define extrapolation state.

◆ extrapolate()

void extrapolate ( int  pdgCode,
double  tof,
const G4ThreeVector &  position,
const G4ThreeVector &  momentum,
const G4ErrorSymMatrix &  covariance 
)

Performs track extrapolation for a single track (specified in genfit2 units).

Parameters
pdgCodeSigned PDG identifier of the particle hypothesis to be used for the extrapolation.
tofStarting time, i.e., time of flight from the IP, at the starting point (ns).
positionStarting point of the extrapolation (cm).
momentumMomentum of the track at the starting point (GeV/c).
covariancePhase-space covariance matrix (6x6) at the starting point (cm, GeV/c).

Definition at line 439 of file TrackExtrapolateG4e.cc.

445{
446 bool isCosmic = false; // DIVOT
447 if ((!m_ExtInitialized) && (!m_MuidInitialized)) {
448 // No EXT nor MUID module in analysis path ==> mimic ext::initialize() with reasonable defaults.
449 // The default values are taken from the EXT module's parameter definitions.
450 Simulation::ExtManager* extMgr = Simulation::ExtManager::GetManager();
451 extMgr->Initialize("Ext", "default", 0.0, 0.25, false, 0, std::vector<std::string>());
452 // Redefine geant4e step length, magnetic field step limitation (fraction of local curvature radius),
453 // and kinetic energy loss limitation (maximum fractional energy loss) by communicating with
454 // the geant4 UI. (Commands were defined in ExtMessenger when physics list was set up.)
455 // *NOTE* If module muid runs after this, its G4UImanager commands will override these.
456 G4UImanager::GetUIpointer()->ApplyCommand("/geant4e/limits/stepLength 250 mm");
457 G4UImanager::GetUIpointer()->ApplyCommand("/geant4e/limits/magField 0.001");
458 G4UImanager::GetUIpointer()->ApplyCommand("/geant4e/limits/energyLoss 0.05");
459 m_DefaultHypotheses = new std::vector<Const::ChargedStable>; // not used
460 initialize(0.1, 0.002, *m_DefaultHypotheses);
461 }
462
463 // Put geant4 in proper state (in case this module is in a separate process)
464 if (G4StateManager::GetStateManager()->GetCurrentState() == G4State_Idle) {
465 G4StateManager::GetStateManager()->SetNewState(G4State_GeomClosed);
466 }
467
468 G4ErrorPropagatorData::GetErrorPropagatorData()->SetTarget(m_TargetExt);
469
470 // Do extrapolation for selected hypothesis (pion, electron, muon, kaon, proton,
471 // deuteron) for the selected track until calorimeter exit.
472
473 G4ThreeVector positionG4e = position * CLHEP::cm; // convert from genfit2 units (cm) to geant4 units (mm)
474 G4ThreeVector momentumG4e = momentum * CLHEP::GeV; // convert from genfit2 units (GeV/c) to geant4 units (MeV/c)
475 if (isCosmic)
476 momentumG4e = -momentumG4e;
477 G4ErrorSymMatrix covarianceG4e(5, 0); // in Geant4e units (GeV/c, cm)
478 fromPhasespaceToG4e(momentum, covariance, covarianceG4e);
479 G4String nameG4e("g4e_" + G4ParticleTable::GetParticleTable()->FindParticle(pdgCode)->GetParticleName());
480 G4ErrorFreeTrajState g4eState(nameG4e, positionG4e, momentumG4e, covarianceG4e);
481 ExtState extState = { nullptr, pdgCode, isCosmic, tof, 0.0, // for EXT and MUID
482 momentumG4e.unit(), 0.0, 0, 0, 0, 0, 0, 0, 0, 0, 0, false // for MUID only
483 };
484 swim(extState, g4eState);
485}
static ExtManager * GetManager()
Get pointer to the instance of this singleton class (create if needed)
Definition: ExtManager.cc:72
void initialize(double minPt, double minKE, std::vector< Const::ChargedStable > &hypotheses)
Initialize for track extrapolation by the EXT module.

◆ findBarrelIntersection()

bool findBarrelIntersection ( ExtState extState,
const G4ThreeVector &  oldPosition,
Intersection intersection 
)
private

Find the intersection point of the track with the crossed BKLM plane.

Definition at line 1408 of file TrackExtrapolateG4e.cc.

1409{
1410 // Be generous: allow outward-moving intersection to be in the dead space between
1411 // largest sensitive-volume Z and m_BarrelHalfLength, not necessarily in a geant4 sensitive volume
1412 if (std::fabs(intersection.position.z() - m_OffsetZ) > m_BarrelHalfLength)
1413 return false;
1414 double phi = intersection.position.phi();
1415 if (phi < 0.0)
1416 phi += TWOPI;
1417 if (phi > TWOPI - PI_8)
1418 phi -= TWOPI;
1419 int sector = (int)((phi + PI_8) / M_PI_4);
1420 int section = intersection.position.z() > m_OffsetZ ?
1423 double oldR = oldPosition * m_BarrelSectorPerp[sector];
1424 double newR = intersection.position * m_BarrelSectorPerp[sector];
1425 for (int layer = extState.firstBarrelLayer; layer <= m_OutermostActiveBarrelLayer; ++layer) {
1426 if (newR < m_BarrelModuleMiddleRadius[section][sector][layer]) break;
1427 if (oldR <= m_BarrelModuleMiddleRadius[section][sector][layer]) {
1428 extState.firstBarrelLayer = layer + 1; // ratchet outward for next call's loop starting value
1429 if (extState.firstBarrelLayer > m_OutermostActiveBarrelLayer) extState.escaped = true;
1430 intersection.inBarrel = true;
1431 intersection.isForward = intersection.position.z() > m_OffsetZ;
1432 intersection.layer = layer;
1433 intersection.sector = sector;
1434 return true;
1435 }
1436 }
1437 return false;
1438}
bool escaped
MUID: flag to indicate that the extrapolated track escaped from the KLM.
int firstBarrelLayer
MUID: outermost barrel layer encountered by the extrapolated track in the prior steps.

◆ findEndcapIntersection()

bool findEndcapIntersection ( ExtState extState,
const G4ThreeVector &  oldPosition,
Intersection intersection 
)
private

Find the intersection point of the track with the crossed EKLM plane.

Definition at line 1440 of file TrackExtrapolateG4e.cc.

1441{
1442 // Be generous: allow intersection to be in the dead space between m_EndcapMinR and innermost
1443 // sensitive-volume radius or between outermost sensitive-volume radius and m_EndcapMaxR,
1444 // not necessarily in a geant4 sensitive volume
1445 if (oldPosition.perp() > m_EndcapMaxR)
1446 return false;
1447 if (intersection.position.perp() < m_EndcapMinR)
1448 return false;
1449 double oldZ = std::fabs(oldPosition.z() - m_OffsetZ);
1450 double newZ = std::fabs(intersection.position.z() - m_OffsetZ);
1451 bool isForward = intersection.position.z() > m_OffsetZ;
1452 int outermostLayer = isForward ? m_OutermostActiveForwardEndcapLayer
1454 for (int layer = extState.firstEndcapLayer; layer <= outermostLayer; ++layer) {
1455 if (newZ < m_EndcapModuleMiddleZ[layer])
1456 break;
1457 if (oldZ <= m_EndcapModuleMiddleZ[layer]) {
1458 extState.firstEndcapLayer = layer + 1; // ratchet outward for next call's loop starting value
1459 if (extState.firstEndcapLayer > outermostLayer)
1460 extState.escaped = true;
1461 intersection.inBarrel = false;
1462 intersection.isForward = isForward;
1463 intersection.layer = layer;
1465 isForward ? 2 : 1, intersection.position) - 1;
1466 return true;
1467 }
1468 }
1469 return false;
1470}
int getSectorByPosition(int section, const HepGeom::Point3D< double > &position) const
Get sector by position.
bool isForward(int cellId)
Check whether the crystal is in forward ECL.
int firstEndcapLayer
MUID: outermost barrel layer encountered by the extrapolated track in the prior steps.

◆ findMatchingBarrelHit()

bool findMatchingBarrelHit ( Intersection intersection,
const Track track 
)
private

Find the matching BKLM 2D hit nearest the intersection point of the track with the crossed BKLM plane.

Definition at line 1472 of file TrackExtrapolateG4e.cc.

1474{
1475 G4ThreeVector extPos0(intersection.position);
1476 double diffBestMagSq = 1.0E60;
1477 int bestHit = -1;
1478 int matchingLayer = intersection.layer + 1;
1479 G4ThreeVector n(m_BarrelSectorPerp[intersection.sector]);
1480 for (int h = 0; h < m_klmHit2ds.getEntries(); ++h) {
1481 KLMHit2d* hit = m_klmHit2ds[h];
1482 if (hit->getSubdetector() != KLMElementNumbers::c_BKLM)
1483 continue;
1484 if (hit->getLayer() != matchingLayer)
1485 continue;
1486 if (hit->isOutOfTime())
1487 continue;
1488 if (std::fabs(hit->getTime() - m_MeanDt) > m_MaxDt)
1489 continue;
1490 G4ThreeVector diff(hit->getPositionX() - intersection.position.x(),
1491 hit->getPositionY() - intersection.position.y(),
1492 hit->getPositionZ() - intersection.position.z());
1493 double dn = diff * n; // in cm
1494 if (std::fabs(dn) < 2.0) {
1495 // Hit and extrapolated point are in the same sector
1496 diff -= n * dn;
1497 if (diff.mag2() < diffBestMagSq) {
1498 diffBestMagSq = diff.mag2();
1499 bestHit = h;
1500 extPos0 = intersection.position;
1501 }
1502 } else {
1503 // Accept a nearby hit in adjacent sector
1504 if (std::fabs(dn) > 50.0)
1505 continue;
1506 int sector = hit->getSector() - 1;
1507 int dSector = abs(intersection.sector - sector);
1508 if ((dSector != +1) && (dSector != m_BarrelNSector - 1))
1509 continue;
1510 // Use the normal vector of the adjacent (hit's) sector
1511 G4ThreeVector nHit(m_BarrelSectorPerp[sector]);
1512 int section = intersection.isForward ?
1515 double dn2 = intersection.position * nHit - m_BarrelModuleMiddleRadius[section][sector][intersection.layer];
1516 dn = diff * nHit + dn2;
1517 if (std::fabs(dn) > 1.0)
1518 continue;
1519 // Project extrapolated track to the hit's plane in the adjacent sector
1520 G4ThreeVector extDir(intersection.momentum.unit());
1521 double extDirA = extDir * nHit;
1522 if (std::fabs(extDirA) < 1.0E-6)
1523 continue;
1524 G4ThreeVector projection = extDir * (dn2 / extDirA);
1525 if (projection.mag() > 15.0)
1526 continue;
1527 diff += projection - nHit * dn;
1528 if (diff.mag2() < diffBestMagSq) {
1529 diffBestMagSq = diff.mag2();
1530 bestHit = h;
1531 extPos0 = intersection.position - projection;
1532 }
1533 }
1534 }
1535
1536 if (bestHit >= 0) {
1537 KLMHit2d* hit = m_klmHit2ds[bestHit];
1538 intersection.isForward = (hit->getSection() == 1);
1539 intersection.sector = hit->getSector() - 1;
1540 intersection.time = hit->getTime();
1541 double localVariance[2] = {m_BarrelScintVariance, m_BarrelScintVariance};
1542 if (hit->inRPC()) {
1543 int nStrips = hit->getPhiStripMax() - hit->getPhiStripMin() + 1;
1544 double dn = nStrips - 1.5;
1545 double factor = std::pow((0.9 + 0.4 * dn * dn), 1.5) * 0.60; // measured-in-Belle resolution
1546 localVariance[0] = m_BarrelPhiStripVariance[intersection.layer] * factor;
1547 nStrips = hit->getZStripMax() - hit->getZStripMin() + 1;
1548 dn = nStrips - 1.5;
1549 factor = std::pow((0.9 + 0.4 * dn * dn), 1.5) * 0.55; // measured-in-Belle resolution
1550 localVariance[1] = m_BarrelZStripVariance[intersection.layer] * factor;
1551 } else {
1552 int nStrips = hit->getPhiStripMax() - hit->getPhiStripMin() + 1;
1553 localVariance[0] *= (nStrips * nStrips); // variance inflated for multi-strip hit
1554 nStrips = hit->getZStripMax() - hit->getZStripMin() + 1;
1555 localVariance[1] *= (nStrips * nStrips); // variance inflated for multi-strip hit
1556 }
1557 G4ThreeVector hitPos(hit->getPositionX(), hit->getPositionY(),
1558 hit->getPositionZ());
1559 adjustIntersection(intersection, localVariance, hitPos, extPos0);
1560 if (intersection.chi2 >= 0.0) {
1561 intersection.hit = bestHit;
1562 hit->isOnTrack(true);
1563 if (track != nullptr) {
1564 track->addRelationTo(hit, intersection.chi2);
1565 RecoTrack* recoTrack = track->getRelatedTo<RecoTrack>();
1567 recoTrack->addBKLMHit(hit, recoTrack->getNumberOfTotalHits() + 1);
1568 }
1569 }
1570 }
1571 }
1572 return intersection.chi2 >= 0.0;
1573
1574}
KLM 2d hit.
Definition: KLMHit2d.h:33
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
bool addBKLMHit(const UsedBKLMHit *bklmHit, const unsigned int sortingParameter, OriginTrackFinder foundByTrackFinder=OriginTrackFinder::c_undefinedTrackFinder)
Adds a bklm hit with the given information to the reco track.
Definition: RecoTrack.h:286
unsigned int getNumberOfTotalHits() const
Return the number of cdc + svd + pxd + bklm + eklm hits.
Definition: RecoTrack.h:436
bool m_addHitsToRecoTrack
Parameter to add the found hits also to the reco tracks or not. Is turned off by default.
void adjustIntersection(Intersection &, const double *, const G4ThreeVector &, const G4ThreeVector &)
Nudge the track using the matching hit.

◆ findMatchingEndcapHit()

bool findMatchingEndcapHit ( Intersection intersection,
const Track track 
)
private

Find the matching EKLM 2D hit nearest the intersection point of the track with the crossed EKLM plane.

Definition at line 1576 of file TrackExtrapolateG4e.cc.

1577{
1578 double diffBestMagSq = 1.0E60;
1579 int bestHit = -1;
1580 int matchingLayer = intersection.layer + 1;
1581 int matchingEndcap = (intersection.isForward ? 2 : 1);
1582 G4ThreeVector n(0.0, 0.0, (intersection.isForward ? 1.0 : -1.0));
1583 for (int h = 0; h < m_klmHit2ds.getEntries(); ++h) {
1584 KLMHit2d* hit = m_klmHit2ds[h];
1585 if (hit->getSubdetector() != KLMElementNumbers::c_EKLM)
1586 continue;
1587 if (hit->getLayer() != matchingLayer)
1588 continue;
1589 if (hit->getSection() != matchingEndcap)
1590 continue;
1591 // DIVOT no such function for EKLM!
1592 // if (hit->isOutOfTime()) continue;
1593 if (std::fabs(hit->getTime() - m_MeanDt) > m_MaxDt)
1594 continue;
1595 G4ThreeVector diff(hit->getPositionX() - intersection.position.x(),
1596 hit->getPositionY() - intersection.position.y(),
1597 hit->getPositionZ() - intersection.position.z());
1598 double dn = diff * n; // in cm
1599 if (std::fabs(dn) > 2.0)
1600 continue;
1601 diff -= n * dn;
1602 if (diff.mag2() < diffBestMagSq) {
1603 diffBestMagSq = diff.mag2();
1604 bestHit = h;
1605 }
1606 }
1607
1608 if (bestHit >= 0) {
1609 KLMHit2d* hit = m_klmHit2ds[bestHit];
1610 intersection.hit = bestHit;
1611 intersection.isForward = (hit->getSection() == 2);
1612 intersection.sector = hit->getSector() - 1;
1613 intersection.time = hit->getTime();
1614 double localVariance[2] = {m_EndcapScintVariance, m_EndcapScintVariance};
1615 int nStrips = hit->getXStripMax() - hit->getXStripMin() + 1;
1616 localVariance[0] *= (nStrips * nStrips); // variance inflated for multi-strip hit
1617 nStrips = hit->getYStripMax() - hit->getYStripMin() + 1;
1618 localVariance[1] *= (nStrips * nStrips); // variance inflated for multi-strip hit
1619 G4ThreeVector hitPos(hit->getPositionX(), hit->getPositionY(), hit->getPositionZ());
1620 adjustIntersection(intersection, localVariance, hitPos, intersection.position);
1621 if (intersection.chi2 >= 0.0) {
1622 // DIVOT no such function for EKLM!
1623 // hit->isOnTrack(true);
1624 if (track != nullptr) {
1625 track->addRelationTo(hit, intersection.chi2);
1626 RecoTrack* recoTrack = track->getRelatedTo<RecoTrack>();
1628 for (const EKLMAlignmentHit& alignmentHit : hit->getRelationsFrom<EKLMAlignmentHit>()) {
1629 recoTrack->addEKLMHit(&alignmentHit, recoTrack->getNumberOfTotalHits() + 1);
1630 }
1631 }
1632 }
1633 }
1634 }
1635 return intersection.chi2 >= 0.0;
1636
1637}
This dataobject is used only for EKLM alignment.
bool addEKLMHit(const UsedEKLMHit *eklmHit, const unsigned int sortingParameter, OriginTrackFinder foundByTrackFinder=OriginTrackFinder::c_undefinedTrackFinder)
Adds an eklm hit with the given information to the reco track.
Definition: RecoTrack.h:300

◆ finishTrack()

void finishTrack ( const ExtState extState,
KLMMuidLikelihood klmMuidLikelihood,
bool  isForward 
)
private

Complete muon identification after end of track extrapolation.

Definition at line 1762 of file TrackExtrapolateG4e.cc.

1763{
1764 /* Done with this track: compute KLM likelihoods and fill the relative dataobject. */
1765 int lastExtLayer = extState.lastBarrelExtLayer + extState.lastEndcapExtLayer + 1;
1767 isForward, extState.escaped, extState.lastBarrelExtLayer, extState.lastEndcapExtLayer);
1768 klmMuidLikelihood->setOutcome(outcome);
1769 klmMuidLikelihood->setIsForward(isForward);
1770 klmMuidLikelihood->setBarrelExtLayer(extState.lastBarrelExtLayer);
1771 klmMuidLikelihood->setEndcapExtLayer(extState.lastEndcapExtLayer);
1772 klmMuidLikelihood->setBarrelHitLayer(extState.lastBarrelHitLayer);
1773 klmMuidLikelihood->setEndcapHitLayer(extState.lastEndcapHitLayer);
1774 klmMuidLikelihood->setExtLayer(lastExtLayer);
1775 klmMuidLikelihood->setHitLayer(((extState.lastEndcapHitLayer == -1) ?
1776 extState.lastBarrelHitLayer :
1777 extState.lastBarrelExtLayer + extState.lastEndcapHitLayer + 1));
1778 klmMuidLikelihood->setChiSquared(extState.chi2);
1779 klmMuidLikelihood->setDegreesOfFreedom(extState.nPoint);
1780 klmMuidLikelihood->setExtLayerPattern(extState.extLayerPattern);
1781 klmMuidLikelihood->setHitLayerPattern(extState.hitLayerPattern);
1782 /* Do KLM likelihood calculation. */
1783 if (outcome != MuidElementNumbers::c_NotReached) { /* Extrapolation reached KLM sensitive volume. */
1784 double denom = 0.0;
1785 int charge = klmMuidLikelihood->getCharge();
1786 std::vector<int> signedPdgVector = MuidElementNumbers::getPDGVector(charge);
1787 std::map<int, double> mapPdgPDF;
1788 for (int pdg : signedPdgVector) {
1789 auto search = m_MuidBuilderMap.find(pdg);
1790 if (search == m_MuidBuilderMap.end())
1791 B2FATAL("Something went wrong: PDF for PDG code " << pdg << " not found!");
1792 double pdf = (search->second)->getPDF(klmMuidLikelihood);
1793 denom += pdf;
1794 mapPdgPDF.insert(std::pair<int, double>(std::abs(pdg), pdf));
1795 }
1796 if (denom < 1.0E-20)
1797 klmMuidLikelihood->setJunkPDFValue(true); /* Anomaly: should be very rare. */
1798 else {
1799 for (auto const& [pdg, pdf] : mapPdgPDF) {
1800 klmMuidLikelihood->setPDFValue(pdf, std::abs(pdg));
1801 if (pdf > 0.0)
1802 klmMuidLikelihood->setLogL(std::log(pdf), std::abs(pdg));
1803 }
1804 }
1805 }
1806}
void setExtLayer(int layer)
Set the outermost EKLM layer crossed in the extrapolation.
void setExtLayerPattern(unsigned int pattern)
Set the pattern of the layers crossed in the extrapolation.
void setLogL(double logL, int pdg)
Set the log-likelihood.
void setChiSquared(double chiSquared)
Set the chi-squared of the extrapolation.
void setHitLayerPattern(unsigned int pattern)
Set the pattern of the layers actually crossed by the track.
void setBarrelExtLayer(int layer)
Set the outermost BKLM layer crossed in the extrapolation.
void setIsForward(bool isForward)
Set if this extrapolation is in forward or backward B/EKLM.
void setJunkPDFValue(bool flag)
Set the junk flag (1 if junk, 0 if not).
void setDegreesOfFreedom(int dof)
Set the number of degrees of freedom (= 2 times the number of KLM hits) for the chi-square computatio...
int getCharge() const
Get the charge of the particle hypothesis used during the extrapolation.
void setPDFValue(double pdfValue, int pdg)
Set the normalized PDF.
void setEndcapExtLayer(int layer)
Set the outermost EKLM layer crossed in the extrapolation.
void setBarrelHitLayer(int layer)
Set the outermost BKLM layer actually crossed by the track.
void setOutcome(unsigned int outcome)
Set the outcome of this extrapolation.
void setHitLayer(int layer)
Set the outermost KLM layer actually crossed by the track.
void setEndcapHitLayer(int layer)
Set the outermost EKLM layer actually crossed by the track.
static unsigned int calculateExtrapolationOutcome(bool isForward, bool escaped, int lastBarrelLayer, int lastEndcapLayer)
Calculate the track extrapolation outcome.
double charge(int pdgCode)
Returns electric charge of a particle with given pdg code.
Definition: EvtPDLUtil.cc:44

◆ fromG4eToPhasespace()

void fromG4eToPhasespace ( const G4ErrorFreeTrajState &  g4eState,
G4ErrorSymMatrix &  covariance 
)
private

Convert the geant4e 5x5 covariance to phasespace 6x6 covariance.

Definition at line 1031 of file TrackExtrapolateG4e.cc.

1032{
1033
1034 // Convert Geant4e covariance matrix with parameters 1/p, lambda, phi, yT, zT (in GeV/c, radians, cm)
1035 // to phase-space covariance matrix with parameters x, y, z, px, py, pz (in GeV/c, cm)
1036 // (1/p) = 1/sqrt( px^2 + py^2 + pz^2 )
1037 // phi = atan( py / px )
1038 // lambda = asin( pz / sqrt( px^2 + py^2 + pz^2 )
1039 // xT = x * cos(lambda) * cos(phi) + y * cos(lambda) * sin(phi) + z * sin(lambda)
1040 // yT = -x * sin(phi) + y * cos(phi)
1041 // zT = -x * sin(lambda) * cos(phi) - y * sin(lambda) * sin(phi) + z * cos(lambda)
1042
1043 G4ErrorFreeTrajParam param = g4eState.GetParameters();
1044 double p = 1.0 / (param.GetInvP() * CLHEP::GeV); // in GeV/c
1045 double pSq = p * p;
1046 double lambda = param.GetLambda(); // in radians
1047 double sinLambda = std::sin(lambda);
1048 double cosLambda = std::cos(lambda);
1049 double phi = param.GetPhi(); // in radians
1050 double sinPhi = std::sin(phi);
1051 double cosPhi = std::cos(phi);
1052
1053 // Transformation Jacobian 6x5 from Geant4e 5x5 to phase-space 6x6
1054
1055 G4ErrorMatrix jacobian(6, 5, 0); // All entries are initialized to 0
1056
1057 jacobian(4, 1) = -pSq * cosLambda * cosPhi; // @(px)/@(1/p)
1058 jacobian(5, 1) = -pSq * cosLambda * sinPhi; // @(py)/@(1/p)
1059 jacobian(6, 1) = -pSq * sinLambda; // @(pz)/@(1/p)
1060
1061 jacobian(4, 2) = -p * sinLambda * cosPhi; // @(px)/@(lambda)
1062 jacobian(5, 2) = -p * sinLambda * sinPhi; // @(py)/@(lambda)
1063 jacobian(6, 2) = p * cosLambda; // @(pz)/@(lambda)
1064
1065 jacobian(4, 3) = -p * cosLambda * sinPhi; // @(px)/@(phi)
1066 jacobian(5, 3) = p * cosLambda * cosPhi; // @(py)/@(phi)
1067
1068 jacobian(1, 4) = -sinPhi; // @(x)/@(yT)
1069 jacobian(2, 4) = cosPhi; // @(y)/@(yT)
1070
1071 jacobian(1, 5) = -sinLambda * cosPhi; // @(x)/@(zT)
1072 jacobian(2, 5) = -sinLambda * sinPhi; // @(y)/@(zT)
1073 jacobian(3, 5) = cosLambda; // @(z)/@(zT)
1074
1075 G4ErrorTrajErr g4eCov = g4eState.GetError();
1076 covariance.assign(g4eCov.similarity(jacobian));
1077
1078}

◆ fromPhasespaceToG4e() [1/2]

void fromPhasespaceToG4e ( const G4ThreeVector &  momentum,
const G4ErrorSymMatrix &  covariance,
G4ErrorTrajErr &  covG4e 
)
private

Convert the phasespace covariance to geant4e covariance.

Definition at line 1133 of file TrackExtrapolateG4e.cc.

1135{
1136
1137 // Convert phase-space covariance matrix with parameters x, y, z, px, py, pz (in genfit2 units: cm, GeV/c)
1138 // to Geant4e covariance matrix with parameters 1/p, lambda, phi, yT, zT (in geant4e units: GeV/c, radians, cm)
1139 // xT = x * cos(lambda) * cos(phi) + y * cos(lambda) * sin(phi) + z * sin(lambda)
1140 // yT = -x * sin(phi) + y * cos(phi)
1141 // zT = -x * sin(lambda) * cos(phi) - y * sin(lambda) * sin(phi) + z * cos(lambda)
1142 // (1/p) = 1/sqrt( px^2 + py^2 + pz^2 )
1143 // phi = atan( py / px )
1144 // lambda = asin( pz / sqrt( px^2 + py^2 + pz^2 )
1145
1146 G4ErrorSymMatrix temp(covariance);
1147
1148 double pInvSq = 1.0 / momentum.mag2();
1149 double pInv = std::sqrt(pInvSq);
1150 double pPerpInv = 1.0 / momentum.perp();
1151 double sinLambda = momentum.cosTheta();
1152 double cosLambda = std::sqrt(1.0 - sinLambda * sinLambda);
1153 double phi = momentum.phi();
1154 double cosPhi = std::cos(phi);
1155 double sinPhi = std::sin(phi);
1156
1157 // Transformation Jacobian 5x6 from phase-space 6x6 to Geant4e 5x5
1158 G4ErrorMatrix jacobian(5, 6, 0);
1159
1160 jacobian(1, 4) = -pInvSq * cosLambda * cosPhi; // @(1/p)/@(px)
1161 jacobian(1, 5) = -pInvSq * cosLambda * sinPhi; // @(1/p)/@(py)
1162 jacobian(1, 6) = -pInvSq * sinLambda; // @(1/p)/@(pz)
1163
1164 jacobian(2, 4) = -pInv * sinLambda * cosPhi; // @(lambda)/@(px)
1165 jacobian(2, 5) = -pInv * sinLambda * sinPhi; // @(lambda)/@(py)
1166 jacobian(2, 6) = pInv * cosLambda; // @(lambda)/@(pz)
1167
1168 jacobian(3, 4) = -pPerpInv * sinPhi; // @(phi)/@(px)
1169 jacobian(3, 5) = pPerpInv * cosPhi; // @(phi)/@(py)
1170
1171 jacobian(4, 1) = -sinPhi; // @(yT)/@(x)
1172 jacobian(4, 2) = cosPhi; // @(yT)/@(y)
1173
1174 jacobian(5, 1) = -sinLambda * cosPhi; // @(zT)/@(x)
1175 jacobian(5, 2) = -sinLambda * sinPhi; // @(zT)/@(y)
1176 jacobian(5, 3) = cosLambda; // @(zT)/@(z)
1177 covG4e = temp.similarity(jacobian);
1178
1179}

◆ fromPhasespaceToG4e() [2/2]

void fromPhasespaceToG4e ( const TVector3 &  momentum,
const TMatrixDSym &  covariance,
G4ErrorTrajErr &  covG4e 
)
private

Convert the phasespace covariance to geant4e covariance.

Definition at line 1080 of file TrackExtrapolateG4e.cc.

1081{
1082
1083 // Convert phase-space covariance matrix with parameters x, y, z, px, py, pz (in genfit2 units: cm, GeV/c)
1084 // to Geant4e covariance matrix with parameters 1/p, lambda, phi, yT, zT (in geant4e units: GeV/c, radians, cm)
1085 // xT = x * cos(lambda) * cos(phi) + y * cos(lambda) * sin(phi) + z * sin(lambda)
1086 // yT = -x * sin(phi) + y * cos(phi)
1087 // zT = -x * sin(lambda) * cos(phi) - y * sin(lambda) * sin(phi) + z * cos(lambda)
1088 // (1/p) = 1/sqrt( px^2 + py^2 + pz^2 )
1089 // phi = atan( py / px )
1090 // lambda = asin( pz / sqrt( px^2 + py^2 + pz^2 )
1091
1092 G4ErrorSymMatrix temp(6, 0);
1093 for (int k = 0; k < 6; ++k) {
1094 for (int j = k; j < 6; ++j) {
1095 temp[j][k] = covariance[j][k];
1096 }
1097 }
1098
1099 double pInvSq = 1.0 / momentum.Mag2();
1100 double pInv = std::sqrt(pInvSq);
1101 double pPerpInv = 1.0 / momentum.Perp();
1102 double sinLambda = momentum.CosTheta();
1103 double cosLambda = std::sqrt(1.0 - sinLambda * sinLambda);
1104 double phi = momentum.Phi();
1105 double cosPhi = std::cos(phi);
1106 double sinPhi = std::sin(phi);
1107
1108 // Transformation Jacobian 5x6 from phase-space 6x6 to Geant4e 5x5
1109 G4ErrorMatrix jacobian(5, 6, 0); // All entries are initialized to 0
1110
1111 jacobian(1, 4) = -pInvSq * cosLambda * cosPhi; // @(1/p)/@(px)
1112 jacobian(1, 5) = -pInvSq * cosLambda * sinPhi; // @(1/p)/@(py)
1113 jacobian(1, 6) = -pInvSq * sinLambda; // @(1/p)/@(pz)
1114
1115 jacobian(2, 4) = -pInv * sinLambda * cosPhi; // @(lambda)/@(px)
1116 jacobian(2, 5) = -pInv * sinLambda * sinPhi; // @(lambda)/@(py)
1117 jacobian(2, 6) = pInv * cosLambda; // @(lambda)/@(pz)
1118
1119 jacobian(3, 4) = -pPerpInv * sinPhi; // @(phi)/@(px)
1120 jacobian(3, 5) = pPerpInv * cosPhi; // @(phi)/@(py)
1121
1122 jacobian(4, 1) = -sinPhi; // @(yT)/@(x)
1123 jacobian(4, 2) = cosPhi; // @(yT)/@(y)
1124
1125 jacobian(5, 1) = -sinLambda * cosPhi; // @(zT)/@(x)
1126 jacobian(5, 2) = -sinLambda * sinPhi; // @(zT)/@(y)
1127 jacobian(5, 3) = cosLambda; // @(zT)/@(z)
1128
1129 covG4e = temp.similarity(jacobian);
1130
1131}

◆ getInstance()

TrackExtrapolateG4e * getInstance ( )
static

Get the singleton's address.

Definition at line 72 of file TrackExtrapolateG4e.cc.

73{
74 if (m_Singleton == nullptr)
76 return m_Singleton;
77}
static TrackExtrapolateG4e * m_Singleton
Stores pointer to the singleton class.
TrackExtrapolateG4e()
constructor is hidden; user calls TrackExtrapolateG4e::getInstance() instead

◆ getStartPoint()

ExtState getStartPoint ( const Track b2track,
int  pdgCode,
G4ErrorFreeTrajState &  g4eState 
)
private

Get the start point for a new reconstructed track with specific PDG hypothesis.

Definition at line 938 of file TrackExtrapolateG4e.cc.

939{
940 ExtState extState = {&b2track, pdgCode, false, 0.0, 0.0, // for EXT and MUID
941 G4ThreeVector(0, 0, 1), 0.0, 0, 0, 0, -1, -1, -1, -1, 0, 0, false // for MUID only
942 };
943 RecoTrack* recoTrack = b2track.getRelatedTo<RecoTrack>();
944 if (recoTrack == nullptr) {
945 B2WARNING("Track without associated RecoTrack: skipping extrapolation for this track.");
946 extState.pdgCode = 0; // prevent start of extrapolation in swim()
947 return extState;
948 }
949 const genfit::AbsTrackRep* trackRep = recoTrack->getCardinalRepresentation();
950 // check for a valid track fit
951 if (!recoTrack->wasFitSuccessful(trackRep)) {
952 B2WARNING("RecoTrack fit failed for cardinal representation: skipping extrapolation for this track.");
953 extState.pdgCode = 0; // prevent start of extrapolation in swim()
954 return extState;
955 }
956 int charge = int(trackRep->getPDGCharge());
957 if (charge != 0) {
958 extState.pdgCode *= charge;
959 } else {
960 charge = 1; // should never happen but persist if it does
961 }
962 TVector3 firstPosition, firstMomentum, lastPosition, lastMomentum; // initialized to zeroes
963 TMatrixDSym firstCov(6), lastCov(6); // initialized to zeroes
964 try {
965 const genfit::MeasuredStateOnPlane& firstState = recoTrack->getMeasuredStateOnPlaneFromFirstHit(trackRep);
966 trackRep->getPosMomCov(firstState, firstPosition, firstMomentum, firstCov);
967 const genfit::MeasuredStateOnPlane& lastState = recoTrack->getMeasuredStateOnPlaneFromLastHit(trackRep);
968 trackRep->getPosMomCov(lastState, lastPosition, lastMomentum, lastCov);
969 // in genfit units (cm, GeV/c)
970 extState.tof = lastState.getTime(); // DIVOT: must be revised when IP profile (reconstructed beam spot) become available!
971 if (lastPosition.Mag2() < firstPosition.Mag2()) {
972 firstPosition = lastPosition;
973 firstMomentum = -lastMomentum;
974 firstCov = lastCov;
975 trackRep->getPosMomCov(firstState, lastPosition, lastMomentum, lastCov);
976 lastMomentum *= -1.0; // extrapolate backwards instead of forwards
977 extState.isCosmic = true;
978 extState.tof = firstState.getTime(); // DIVOT: must be revised when IP profile (reconstructed beam spot) become available!
979 }
980
981 G4ParticleDefinition* particle = G4ParticleTable::GetParticleTable()->FindParticle(extState.pdgCode);
982 if (extState.pdgCode != trackRep->getPDG()) {
983 double pSq = lastMomentum.Mag2();
984 double mass = particle->GetPDGMass() / CLHEP::GeV;
985 extState.tof *= std::sqrt((pSq + mass * mass) / (pSq + lastState.getMass() * lastState.getMass()));
986 }
987
988 extState.directionAtIP.set(firstMomentum.Unit().X(), firstMomentum.Unit().Y(), firstMomentum.Unit().Z());
989 if (m_MagneticField != 0.0) { // in gauss
990 double radius = (firstMomentum.Perp() * CLHEP::GeV / CLHEP::eV) / (CLHEP::c_light * charge * m_MagneticField); // in cm
991 double centerPhi = extState.directionAtIP.phi() - M_PI_2;
992 double centerX = firstPosition.X() + radius * std::cos(centerPhi);
993 double centerY = firstPosition.Y() + radius * std::sin(centerPhi);
994 double pocaPhi = atan2(charge * centerY, charge * centerX) + M_PI;
995 double ipPerp = extState.directionAtIP.perp();
996 if (ipPerp > 0.0) {
997 extState.directionAtIP.setX(+std::sin(pocaPhi) * ipPerp);
998 extState.directionAtIP.setY(-std::cos(pocaPhi) * ipPerp);
999 }
1000 } else {
1001 // No field: replace flaky covariance matrix with a diagonal one measured in 1.5T field
1002 // for a 10 GeV/c track ... and replace momentum magnitude with fixed 10 GeV/c
1003 lastCov *= 0.0;
1004 lastCov[0][0] = 5.0E-5;
1005 lastCov[1][1] = 1.0E-7;
1006 lastCov[2][2] = 5.0E-4;
1007 lastCov[3][3] = 3.5E-3;
1008 lastCov[4][4] = 3.5E-3;
1009 lastMomentum = lastMomentum.Unit() * 10.0;
1010 }
1011
1012 G4ThreeVector posG4e(lastPosition.X() * CLHEP::cm, lastPosition.Y() * CLHEP::cm,
1013 lastPosition.Z() * CLHEP::cm); // in Geant4 units (mm)
1014 G4ThreeVector momG4e(lastMomentum.X() * CLHEP::GeV, lastMomentum.Y() * CLHEP::GeV,
1015 lastMomentum.Z() * CLHEP::GeV); // in Geant4 units (MeV/c)
1016 G4ErrorSymMatrix covG4e(5, 0); // in Geant4e units (GeV/c, cm)
1017 fromPhasespaceToG4e(lastMomentum, lastCov, covG4e); // in Geant4e units (GeV/c, cm)
1018 g4eState.SetData("g4e_" + particle->GetParticleName(), posG4e, momG4e);
1019 g4eState.SetParameters(posG4e, momG4e); // compute private-state parameters from momG4e
1020 g4eState.SetError(covG4e);
1021 } catch (const genfit::Exception&) {
1022 B2WARNING("genfit::MeasuredStateOnPlane() exception: skipping extrapolation for this track. initial momentum = ("
1023 << firstMomentum.X() << "," << firstMomentum.Y() << "," << firstMomentum.Z() << ")");
1024 extState.pdgCode = 0; // prevent start of extrapolation in swim()
1025 }
1026
1027 return extState;
1028}
bool wasFitSuccessful(const genfit::AbsTrackRep *representation=nullptr) const
Returns true if the last fit with the given representation was successful.
Definition: RecoTrack.cc:336
genfit::AbsTrackRep * getCardinalRepresentation() const
Get a pointer to the cardinal track representation. You are not allowed to modify or delete it!
Definition: RecoTrack.h:631
const genfit::MeasuredStateOnPlane & getMeasuredStateOnPlaneFromLastHit(const genfit::AbsTrackRep *representation=nullptr) const
Return genfit's MeasuredStateOnPlane for the last hit in a fit useful for extrapolation of measuremen...
Definition: RecoTrack.cc:619
const genfit::MeasuredStateOnPlane & getMeasuredStateOnPlaneFromFirstHit(const genfit::AbsTrackRep *representation=nullptr) const
Return genfit's MeasuredStateOnPlane for the first hit in a fit useful for extrapolation of measureme...
Definition: RecoTrack.cc:605
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
G4ThreeVector directionAtIP
MUID: initial direction of track, used for KLID.

◆ getVolumeID()

void getVolumeID ( const G4TouchableHandle &  touch,
Const::EDetector detID,
int &  copyID 
)
private

Get the physical volume information for a geant4 physical volume.

Definition at line 852 of file TrackExtrapolateG4e.cc.

853{
854
855 // default values
856 detID = Const::EDetector::invalidDetector;
857 copyID = 0;
858
859 G4VPhysicalVolume* pv = touch->GetVolume(0);
860 std::map<G4VPhysicalVolume*, enum VolTypes>::iterator it = m_EnterExit->find(pv);
861 if (it == m_EnterExit->end())
862 return;
863
864 switch (it->second) {
865 case VOLTYPE_CDC:
866 detID = Const::EDetector::CDC;
867 copyID = pv->GetCopyNo();
868 return;
869 case VOLTYPE_TOP1:
870 detID = Const::EDetector::TOP;
871 copyID = -(pv->GetCopyNo()); // negative to distinguish module and quartz hits
872 return;
873 case VOLTYPE_TOP2:
874 detID = Const::EDetector::TOP;
875 if (touch->GetHistoryDepth() >= 1)
876 copyID = touch->GetVolume(1)->GetCopyNo();
877 return;
878 case VOLTYPE_TOP3:
879 detID = Const::EDetector::TOP;
880 if (touch->GetHistoryDepth() >= 2)
881 copyID = touch->GetVolume(2)->GetCopyNo();
882 return;
883 case VOLTYPE_ARICH1:
884 detID = Const::EDetector::ARICH;
885 copyID = 12345;
886 return;
887 case VOLTYPE_ARICH2:
888 detID = Const::EDetector::ARICH;
889 copyID = 6789;
890 return;
891 case VOLTYPE_ARICH3:
892 detID = Const::EDetector::ARICH;
893 if (touch->GetHistoryDepth() >= 2)
894 copyID = touch->GetVolume(2)->GetCopyNo();
895 return;
896 case VOLTYPE_ECL:
897 detID = Const::EDetector::ECL;
899 return;
900 case VOLTYPE_BKLM1: // BKLM RPCs
901 detID = Const::EDetector::BKLM;
902 if (touch->GetHistoryDepth() == DEPTH_RPC) {
903 // int plane = touch->GetCopyNumber(0);
904 int layer = touch->GetCopyNumber(4);
905 int sector = touch->GetCopyNumber(6);
906 int section = touch->GetCopyNumber(7);
907 copyID = BKLMElementNumbers::moduleNumber(section, sector, layer);
908 }
909 return;
910 case VOLTYPE_BKLM2: // BKLM scints
911 detID = Const::EDetector::BKLM;
912 if (touch->GetHistoryDepth() == DEPTH_SCINT) {
913 int strip = touch->GetCopyNumber(1);
914 int plane = (touch->GetCopyNumber(2) == BKLM_INNER) ?
917 int layer = touch->GetCopyNumber(6);
918 int sector = touch->GetCopyNumber(8);
919 int section = touch->GetCopyNumber(9);
921 section, sector, layer, plane, strip);
922 BKLMStatus::setMaximalStrip(copyID, strip);
923 }
924 return;
925 case VOLTYPE_EKLM:
926 detID = Const::EDetector::EKLM;
928 touch->GetVolume(7)->GetCopyNo(),
929 touch->GetVolume(6)->GetCopyNo(),
930 touch->GetVolume(5)->GetCopyNo(),
931 touch->GetVolume(4)->GetCopyNo(),
932 touch->GetVolume(1)->GetCopyNo());
933 return;
934 }
935
936}
static KLMChannelNumber channelNumber(int section, int sector, int layer, int plane, int strip)
Get channel number.
static KLMModuleNumber moduleNumber(int section, int sector, int layer, bool fatalError=true)
Get module number.
static void setMaximalStrip(int &module, int strip)
Set maximal strip number.
Definition: BKLMStatus.h:49
static ECLGeometryPar * Instance()
Static method to get a reference to the ECLGeometryPar instance.
int ECLVolumeToCellID(const G4VTouchable *)
Get Cell Id (LEP: new way)
int stripNumber(int section, int layer, int sector, int plane, int strip) const
Get strip number.
@ VOLTYPE_ARICH2
ARICH Img plate.
@ VOLTYPE_CDC
CDC.
@ VOLTYPE_TOP2
TOP quartz.
@ VOLTYPE_TOP3
TOP glue.
@ VOLTYPE_ARICH3
ARICH HAPD window.
@ VOLTYPE_BKLM2
BKLM scintillator.
@ VOLTYPE_EKLM
EKLM.
@ VOLTYPE_BKLM1
BKLM RPC.
@ VOLTYPE_ARICH1
ARICH aerogel.
@ VOLTYPE_ECL
ECL.
@ VOLTYPE_TOP1
TOP container.

◆ initialize() [1/2]

void initialize ( double  meanDt,
double  maxDt,
double  maxSeparation,
double  maxKLMTrackClusterDistance,
double  maxECLTrackClusterDistance,
double  minPt,
double  minKE,
bool  addHitsToRecoTrack,
std::vector< Const::ChargedStable > &  hypotheses 
)

Initialize for track extrapolation by the MUID module.

Parameters
meanDtMean value of the in-time window (ns).
maxDtHalf-width of the in-time window (ns).
maxSeparationMaximum separation between track crossing and matching hit in detector plane (number of sigmas).
maxKLMTrackClusterDistanceMaximum distance between associated track and KLMCluster (cm), criterion for matching relation Track->KLMCluster on MDST.
maxECLTrackClusterDistanceMaximum distance between associated track and ECLCluster (cm).
minPtMinimum transverse momentum to begin extrapolation (GeV/c).
minKEMinimum kinetic energy to continue extrapolation (GeV/c).
addHitsToRecoTrackParameter to add the found hits also to the reco tracks or not. Is turned off by default.
hypothesesVector of charged-particle hypotheses used in extrapolation of each track.

Definition at line 191 of file TrackExtrapolateG4e.cc.

195{
196 m_MuidInitialized = true;
197 m_addHitsToRecoTrack = addHitsToRecoTrack;
198
199 // Define required objects, register the new ones and relations
200 m_eclClusters.isRequired();
201 m_klmHit2ds.isRequired();
202 m_klmClusters.isRequired();
204 m_tracks.isRequired();
205 m_extHits.registerInDataStore();
206 m_klmMuidLikelihoods.registerInDataStore();
207 m_klmMuidHits.registerInDataStore();
208 m_trackClusterSeparations.registerInDataStore();
209 m_tracks.registerRelationTo(m_extHits);
210 m_tracks.registerRelationTo(m_klmMuidLikelihoods);
211 m_tracks.registerRelationTo(m_klmMuidHits);
212 m_tracks.registerRelationTo(m_klmHit2ds);
213 m_tracks.registerRelationTo(m_trackClusterSeparations);
214 m_tracks.registerRelationTo(m_klmClusters);
215 m_klmClusters.registerRelationTo(m_trackClusterSeparations);
216 m_eclClusters.registerRelationTo(m_extHits);
218
219 // Save the in-time cut's central value and width for valid hits
220 m_MeanDt = meanDt;
221 m_MaxDt = maxDt;
222
223 // Save the magnetic field z component (gauss) at the origin
224 m_MagneticField = BFieldManager::getField(0, 0, 0).Z() / Unit::T * CLHEP::tesla / CLHEP::gauss;
225
226 // Convert from sigma to variance for hit-position uncertainty
227 m_MaxDistSqInVariances = maxKLMTrackHitDistance * maxKLMTrackHitDistance;
228
229 // Convert user cutoff values to geant4 units
230 m_MaxKLMTrackClusterDistance = std::max(0.0, maxKLMTrackClusterDistance) * CLHEP::cm; // mm
231 m_MaxECLTrackClusterDistance = std::max(0.0, maxECLTrackClusterDistance) * CLHEP::cm; // mm
232 m_MinPt = std::max(0.0, minPt) * CLHEP::GeV;
233 m_MinKE = std::max(0.0, minKE) * CLHEP::GeV;
234
235 // Save pointer to the list of particle hypotheses for EXT extrapolation
236 m_HypothesesMuid = &hypotheses;
237
238 // Define the list of volumes that will have their entry and/or
239 // exit points stored during the extrapolation.
241
242 // Store the address of the ExtManager (used later)
244
245 // Set up the EXT-specific geometry (might have already been done by EXT)
246 if (m_TargetExt == nullptr) {
247 if (!m_COILGeometryPar.isValid())
248 B2FATAL("Coil geometry data are not available.");
249 double offsetZ = m_COILGeometryPar->getGlobalOffsetZ();
250 double rMinCoil = m_COILGeometryPar->getCryoRmin();
251 double halfLength = m_COILGeometryPar->getCryoLength();
252 m_COILGeometryPar.addCallback([this, &offsetZ, &rMinCoil, &halfLength]() {
253 offsetZ = m_COILGeometryPar->getGlobalOffsetZ();
254 rMinCoil = m_COILGeometryPar->getCryoRmin();
255 halfLength = m_COILGeometryPar->getCryoLength();
256 });
257 m_TargetExt = new Simulation::ExtCylSurfaceTarget(rMinCoil, offsetZ - halfLength, offsetZ + halfLength);
258 G4ErrorPropagatorData::GetErrorPropagatorData()->SetTarget(m_TargetExt);
259 }
260 if (!m_BeamPipeGeo.isValid())
261 B2FATAL("Beam pipe geometry data are not available.");
262 double beampipeRadius = m_BeamPipeGeo->getParameter("Lv2OutBe.R2") * CLHEP::cm; // mm
263 m_BeamPipeGeo.addCallback([this, &beampipeRadius]() {
264 beampipeRadius = m_BeamPipeGeo->getParameter("Lv2OutBe.R2") * CLHEP::cm;
265 });
266 m_MinRadiusSq = beampipeRadius * beampipeRadius; // mm^2
267
268 // Set up the MUID-specific geometry
270 const EKLM::GeometryData& eklmGeometry = EKLM::GeometryData::Instance();
271 m_BarrelHalfLength = bklmGeometry->getHalfLength() * CLHEP::cm; // in G4 units (mm)
272 m_EndcapHalfLength = 0.5 * eklmGeometry.getSectionPosition()->getLength(); // in G4 units (mm)
273 m_OffsetZ = bklmGeometry->getOffsetZ() * CLHEP::cm; // in G4 units (mm)
274 double minZ = m_OffsetZ - (m_BarrelHalfLength + 2.0 * m_EndcapHalfLength);
275 double maxZ = m_OffsetZ + (m_BarrelHalfLength + 2.0 * m_EndcapHalfLength);
276 m_BarrelNSector = bklmGeometry->getNSector();
277 m_BarrelMaxR = bklmGeometry->getOuterRadius() * CLHEP::cm / std::cos(M_PI / m_BarrelNSector); // in G4 units (mm)
278 m_TargetMuid = new Simulation::ExtCylSurfaceTarget(m_BarrelMaxR, minZ, maxZ);
279 G4ErrorPropagatorData::GetErrorPropagatorData()->SetTarget(m_TargetMuid);
280
281 m_BarrelHalfLength /= CLHEP::cm; // now in G4e units (cm)
282 m_EndcapHalfLength /= CLHEP::cm; // now in G4e units (cm)
283 m_OffsetZ /= CLHEP::cm; // now in G4e units (cm)
284 m_BarrelMinR = bklmGeometry->getGap1InnerRadius(); // in G4e units (cm)
285 m_BarrelMaxR /= CLHEP::cm; // now in G4e units (cm)
286 m_EndcapMinR = eklmGeometry.getSectionPosition()->getInnerR() / CLHEP::cm; // in G4e units (cm)
287 m_EndcapMaxR = eklmGeometry.getSectionPosition()->getOuterR() / CLHEP::cm; // in G4e units (cm)
289
290 // Measurement uncertainties and acceptance windows
291 double width = eklmGeometry.getStripGeometry()->getWidth() / CLHEP::cm; // in G4e units (cm)
292 m_EndcapScintVariance = width * width / 12.0;
293 width = bklmGeometry->getScintHalfWidth() * 2.0; // in G4e units (cm)
294 m_BarrelScintVariance = width * width / 12.0;
295 for (int layer = 1; layer <= BKLMElementNumbers::getMaximalLayerNumber(); ++layer) {
296 const bklm::Module* module =
297 bklmGeometry->findModule(BKLMElementNumbers::c_ForwardSection, 1, layer);
298 width = module->getPhiStripWidth(); // in G4e units (cm)
299 m_BarrelPhiStripVariance[layer - 1] = width * width / 12.0;
300 width = module->getZStripWidth(); // in G4e units (cm)
301 m_BarrelZStripVariance[layer - 1] = width * width / 12.0;
302 }
303
304 // KLM geometry (for associating KLM hit with extrapolated crossing point)
306 for (int sector = 1; sector <= BKLMElementNumbers::getMaximalSectorNumber(); ++sector) {
307 double phi = M_PI_4 * (sector - 1);
308 m_BarrelSectorPerp[sector - 1].set(std::cos(phi), std::sin(phi), 0.0);
309 m_BarrelSectorPhi[sector - 1].set(-std::sin(phi), std::cos(phi), 0.0);
310 }
312 for (KLMChannelIndex& klmLayer : klmLayers) {
313 if (klmLayer.getSubdetector() == KLMElementNumbers::c_BKLM)
314 m_BarrelModuleMiddleRadius[klmLayer.getSection()][klmLayer.getSector() - 1][klmLayer.getLayer() - 1] =
315 bklmGeometry->getActiveMiddleRadius(klmLayer.getSection(), klmLayer.getSector(), klmLayer.getLayer()); // in G4e units (cm)
316 }
317 double dz(eklmGeometry.getLayerShiftZ() / CLHEP::cm); // in G4e units (cm)
318 double z0((eklmGeometry.getSectionPosition()->getZ()
319 + eklmGeometry.getLayerShiftZ()
320 - 0.5 * eklmGeometry.getSectionPosition()->getLength()
321 - 0.5 * eklmGeometry.getLayerPosition()->getLength()) / CLHEP::cm); // in G4e units (cm)
323 1; // zero-based counting
325 1; // zero-based counting
326 for (int layer = 1; layer <= EKLMElementNumbers::getMaximalLayerNumber(); ++layer) {
327 m_EndcapModuleMiddleZ[layer - 1] = z0 + dz * (layer - 1); // in G4e units (cm)
328 }
329
331}
int getMaximalDetectorLayerNumber(int section) const
Get maximal detector layer number.
static constexpr int getMaximalLayerNumber()
Get maximal layer number.
double getOuterR() const
Get outer radius.
Definition: EKLMGeometry.h:129
double getInnerR() const
Get inner radius.
Definition: EKLMGeometry.h:112
double getZ() const
Get Z coordinate.
Definition: EKLMGeometry.h:197
double getLength() const
Get length.
Definition: EKLMGeometry.h:146
double getWidth() const
Get width.
const ElementPosition * getLayerPosition() const
Get position data for layers.
double getLayerShiftZ() const
Get Z distance between two layers.
const StripGeometry * getStripGeometry() const
Get strip geometry data.
const ElementPosition * getSectionPosition() const
Get position data for sections.
EKLM geometry data.
Definition: GeometryData.h:38
static const GeometryData & Instance(enum DataSource dataSource=c_Database, const GearDir *gearDir=nullptr)
Instantiation.
Definition: GeometryData.cc:33
static const TransformDataGlobalAligned & Instance()
Instantiation.
KLM channel index.
static void registerRequiredRelations(StoreArray< RecoTrack > &recoTracks, std::string const &pxdHitsStoreArrayName="", std::string const &svdHitsStoreArrayName="", std::string const &cdcHitsStoreArrayName="", std::string const &bklmHitsStoreArrayName="", std::string const &eklmHitsStoreArrayName="", std::string const &recoHitInformationStoreArrayName="")
Convenience method which registers all relations required to fully use a RecoTrack.
Definition: RecoTrack.cc:53
bool isRequired(const std::string &name="")
Ensure this array/object has been registered previously.
StoreArray< KLMMuidLikelihood > m_klmMuidLikelihoods
KLM muid likelihoods.
StoreArray< TrackClusterSeparation > m_trackClusterSeparations
Track cluster separation.
StoreArray< RecoTrack > m_recoTracks
Reco tracks.
void registerVolumes()
Register the list of geant4 physical volumes whose entry/exit points will be saved during extrapolati...
DBObjPtr< COILGeometryPar > m_COILGeometryPar
Conditions-database object for COIL geometry.
DBObjPtr< BeamPipeGeo > m_BeamPipeGeo
Conditions-database object for beam pipe geometry.
static const double T
[tesla]
Definition: Unit.h:120
Provides BKLM geometry parameters for simulation, reconstruction etc (from Gearbox or DataBase)
Definition: GeometryPar.h:37
double getScintHalfWidth(void) const
Get the height of the entire volume of a scintillator strip (including TiO2 coating)
Definition: GeometryPar.h:118
double getGap1InnerRadius(void) const
Get the radius of the inner tangent circle of gap 0 (innermost)
Definition: GeometryPar.h:211
double getOuterRadius(void) const
Get the radius of the inscribed circle of the outer polygon.
Definition: GeometryPar.h:187
double getOffsetZ(void) const
Get the global shift along a of the entire BKLM.
Definition: GeometryPar.h:157
int getNSector(void) const
Get the number of sectors of the BKLM.
Definition: GeometryPar.h:175
double getHalfLength(void) const
Get the half-length along z of the BKLM.
Definition: GeometryPar.h:181
double getActiveMiddleRadius(int section, int sector, int layer) const
Get the radial midpoint of the detector module's active volume of specified layer.
Definition: GeometryPar.cc:607
static void getField(const double *pos, double *field)
return the magnetic field at a given position.
Definition: BFieldManager.h:91

◆ initialize() [2/2]

void initialize ( double  minPt,
double  minKE,
std::vector< Const::ChargedStable > &  hypotheses 
)

Initialize for track extrapolation by the EXT module.

Parameters
minPtMinimum transverse momentum to begin extrapolation (GeV/c).
minKEMinimum kinetic energy to continue extrapolation (GeV/c).
hypothesesVector of charged-particle hypotheses used in extrapolation of each track.

Definition at line 138 of file TrackExtrapolateG4e.cc.

140{
141 m_ExtInitialized = true;
142
143 // Define required objects, register the new ones and relations
145 m_tracks.isRequired();
146 m_extHits.registerInDataStore();
147 m_tracks.registerRelationTo(m_extHits);
148
149 // Save the magnetic field z component (gauss) at the origin
150 m_MagneticField = BFieldManager::getField(0, 0, 0).Z() / Unit::T * CLHEP::tesla / CLHEP::gauss;
151
152 // Convert user cutoff values to geant4 units
153 m_MinPt = std::max(0.0, minPt) * CLHEP::GeV;
154 m_MinKE = std::max(0.0, minKE) * CLHEP::GeV;
155
156 // Save pointer to the list of particle hypotheses for EXT extrapolation
157 m_HypothesesExt = &hypotheses;
158
159 // Define the list of volumes that will have their entry and/or
160 // exit points stored during the extrapolation.
162
163 // Store the address of the ExtManager (used later)
165
166 // Set up the EXT-specific geometry (might have already been done by MUID)
167 if (m_TargetExt == nullptr) {
168 if (!m_COILGeometryPar.isValid())
169 B2FATAL("Coil geometry data are not available.");
170 double offsetZ = m_COILGeometryPar->getGlobalOffsetZ();
171 double rMinCoil = m_COILGeometryPar->getCryoRmin();
172 double halfLength = m_COILGeometryPar->getCryoLength();
173 m_COILGeometryPar.addCallback([this, &offsetZ, &rMinCoil, &halfLength]() {
174 offsetZ = m_COILGeometryPar->getGlobalOffsetZ();
175 rMinCoil = m_COILGeometryPar->getCryoRmin();
176 halfLength = m_COILGeometryPar->getCryoLength();
177 });
178 m_TargetExt = new Simulation::ExtCylSurfaceTarget(rMinCoil, offsetZ - halfLength, offsetZ + halfLength);
179 G4ErrorPropagatorData::GetErrorPropagatorData()->SetTarget(m_TargetExt);
180 }
181 if (!m_BeamPipeGeo.isValid())
182 B2FATAL("Beam pipe geometry data are not available.");
183 double beampipeRadius = m_BeamPipeGeo->getParameter("Lv2OutBe.R2") * CLHEP::cm; // mm
184 m_BeamPipeGeo.addCallback([this, &beampipeRadius]() {
185 beampipeRadius = m_BeamPipeGeo->getParameter("Lv2OutBe.R2") * CLHEP::cm;
186 });
187 m_MinRadiusSq = beampipeRadius * beampipeRadius; // mm^2
188}

◆ registerVolumes()

void registerVolumes ( )
private

Register the list of geant4 physical volumes whose entry/exit points will be saved during extrapolation.

Definition at line 788 of file TrackExtrapolateG4e.cc.

789{
790 G4PhysicalVolumeStore* pvStore = G4PhysicalVolumeStore::GetInstance();
791 if (pvStore->size() == 0) {
792 B2FATAL("No geometry defined. Please create the geometry first.");
793 }
794 if (m_EnterExit != nullptr) // Only do this once
795 return;
796 m_EnterExit = new std::map<G4VPhysicalVolume*, enum VolTypes>;
797 m_BKLMVolumes = new std::vector<G4VPhysicalVolume*>;
798 for (std::vector<G4VPhysicalVolume*>::iterator iVol = pvStore->begin();
799 iVol != pvStore->end(); ++iVol) {
800 const G4String name = (*iVol)->GetName();
801 // Do not store ExtHits in CDC, so let's start from TOP and ARICH.
802 // TOP doesn't have one envelope; it has 16 "TOPModule"s
803 if (name.find("TOPModule") != std::string::npos) {
804 (*m_EnterExit)[*iVol] = VOLTYPE_TOP1;
805 }
806 // TOP quartz bar (=sensitive)
807 else if (name.find("_TOPPrism_") != std::string::npos ||
808 name.find("_TOPBarSegment") != std::string::npos ||
809 name.find("_TOPMirrorSegment") != std::string::npos) {
810 (*m_EnterExit)[*iVol] = VOLTYPE_TOP2;
811 }
812 // TOP quartz glue (not sensitive?)
813 else if (name.find("TOPBarSegment1Glue") != std::string::npos ||
814 name.find("TOPBarSegment2Glue") != std::string::npos ||
815 name.find("TOPMirrorSegmentGlue") != std::string::npos) {
816 (*m_EnterExit)[*iVol] = VOLTYPE_TOP3;
817 // ARICH volumes
818 } else if (name == "ARICH.AerogelSupportPlate") {
819 (*m_EnterExit)[*iVol] = VOLTYPE_ARICH1;
820 } else if (name == "ARICH.AerogelImgPlate") {
821 (*m_EnterExit)[*iVol] = VOLTYPE_ARICH2;
822 } else if (name.find("ARICH.HAPDWindow") != std::string::npos) {
823 (*m_EnterExit)[*iVol] = VOLTYPE_ARICH3;
824 }
825 // ECL crystal
826 else if (name.find("lv_barrel_crystal_") != std::string::npos ||
827 name.find("lv_forward_crystal_") != std::string::npos ||
828 name.find("lv_backward_crystal_") != std::string::npos) {
829 (*m_EnterExit)[*iVol] = VOLTYPE_ECL;
830 }
831 // Barrel KLM: BKLM.Layer**GasPhysical for RPCs or BKLM.Layer**ChimneyGasPhysical for RPCs
832 // BKLM.ScintActiveType*Physical for scintillator strips
833 else if (name.compare(0, 5, "BKLM.") == 0) {
834 if (name.find("GasPhysical") != std::string::npos) {
835 (*m_EnterExit)[*iVol] = VOLTYPE_BKLM1;
836 } else if (name.find("ScintActiveType") != std::string::npos) {
837 (*m_EnterExit)[*iVol] = VOLTYPE_BKLM2;
838 } else if ((name.find("ScintType") != std::string::npos) ||
839 (name.find("ElectrodePhysical") != std::string::npos)) {
840 m_BKLMVolumes->push_back(*iVol);
841 }
842 }
843 // Endcap KLM: StripSensitive_*
844 else if (name.compare(0, 14, "StripSensitive") == 0) {
845 (*m_EnterExit)[*iVol] = VOLTYPE_EKLM;
846 }
847 }
848
849}

◆ swim() [1/2]

void swim ( ExtState extState,
G4ErrorFreeTrajState &  g4eState 
)
private

Swim a single track (EXT) until it stops or leaves the target cylinder.

Definition at line 707 of file TrackExtrapolateG4e.cc.

708{
709 if (extState.pdgCode == 0)
710 return;
711 if (g4eState.GetMomentum().perp() <= m_MinPt)
712 return;
713 if (m_TargetExt->GetDistanceFromPoint(g4eState.GetPosition()) < 0.0)
714 return;
715 G4ParticleDefinition* particle = G4ParticleTable::GetParticleTable()->FindParticle(extState.pdgCode);
716 double mass = particle->GetPDGMass();
717 double minPSq = (mass + m_MinKE) * (mass + m_MinKE) - mass * mass;
718 G4ErrorMode propagationMode = (extState.isCosmic ? G4ErrorMode_PropBackwards : G4ErrorMode_PropForwards);
719 m_ExtMgr->InitTrackPropagation(propagationMode);
720 while (true) {
721 const G4int errCode = m_ExtMgr->PropagateOneStep(&g4eState, propagationMode);
722 G4Track* track = g4eState.GetG4Track();
723 const G4Step* step = track->GetStep();
724 const G4StepPoint* preStepPoint = step->GetPreStepPoint();
725 const G4StepPoint* postStepPoint = step->GetPostStepPoint();
726 G4TouchableHandle preTouch = preStepPoint->GetTouchableHandle();
727 G4VPhysicalVolume* pVol = preTouch->GetVolume();
728 const G4int preStatus = preStepPoint->GetStepStatus();
729 const G4int postStatus = postStepPoint->GetStepStatus();
730 G4ThreeVector pos = track->GetPosition(); // this is at postStepPoint
731 G4ThreeVector mom = track->GetMomentum(); // ditto
732 // First step on this track?
733 if (extState.isCosmic)
734 mom = -mom;
735 if (preStatus == fUndefined) {
736 if (m_EnterExit->find(pVol) != m_EnterExit->end()) {
737 createExtHit(EXT_FIRST, extState, g4eState, preStepPoint, preTouch);
738 }
739 }
740 // Ignore the zero-length step by PropagateOneStep() at each boundary
741 if (step->GetStepLength() > 0.0) {
742 double dt = step->GetDeltaTime();
743 double dl = step->GetStepLength() / track->GetMaterial()->GetRadlen();
744 if (preStatus == fGeomBoundary) { // first step in this volume?
745 if (m_EnterExit->find(pVol) != m_EnterExit->end()) {
746 createExtHit(EXT_ENTER, extState, g4eState, preStepPoint, preTouch);
747 }
748 }
749 if (extState.isCosmic) {
750 extState.tof -= dt;
751 extState.length -= dl;
752 } else {
753 extState.tof += dt;
754 extState.length += dl;
755 }
756 // Last step in this volume?
757 if (postStatus == fGeomBoundary) {
758 if (m_EnterExit->find(pVol) != m_EnterExit->end()) {
759 createExtHit(EXT_EXIT, extState, g4eState, postStepPoint, preTouch);
760 }
761 }
762 }
763 // Post-step momentum too low?
764 if (errCode || (mom.mag2() < minPSq)) {
765 if (m_EnterExit->find(pVol) != m_EnterExit->end()) {
766 createExtHit(EXT_STOP, extState, g4eState, postStepPoint, preTouch);
767 }
768 break;
769 }
770 // Detect escapes from the imaginary target cylinder.
771 if (m_TargetExt->GetDistanceFromPoint(pos) < 0.0) {
772 if (m_EnterExit->find(pVol) != m_EnterExit->end()) {
773 createExtHit(EXT_ESCAPE, extState, g4eState, postStepPoint, preTouch);
774 }
775 break;
776 }
777 // Stop extrapolating as soon as the track curls inward too much
778 if (pos.perp2() < m_MinRadiusSq) {
779 break;
780 }
781 } // track-extrapolation "infinite" loop
782
783 m_ExtMgr->EventTermination(propagationMode);
784
785}
virtual G4double GetDistanceFromPoint(const G4ThreeVector &point, const G4ThreeVector &direc) const
Get the distance from a point to the cylinder along direc.
void EventTermination(G4ErrorMode)
Terminate an event and set state to G4ErrorState_Init.
Definition: ExtManager.cc:136
void InitTrackPropagation(G4ErrorMode)
Initialize for propagation of a track and set state to G4ErrorState_Propagating.
Definition: ExtManager.cc:106
G4int PropagateOneStep(G4ErrorTrajState *currentTS, G4ErrorMode mode=G4ErrorMode_PropForwards)
Propagate a track by one step.
Definition: ExtManager.cc:122
void createExtHit(ExtHitStatus, const ExtState &, const G4ErrorFreeTrajState &, const G4StepPoint *, const G4TouchableHandle &)
Create another EXT extrapolation hit for a track candidate.
double length
Length from start of extrapolation (rad lengths), updated during extrapolation.

◆ swim() [2/2]

void swim ( ExtState extState,
G4ErrorFreeTrajState &  g4eState,
const std::vector< std::pair< ECLCluster *, G4ThreeVector > > *  eclClusterInfo,
const std::vector< std::pair< KLMCluster *, G4ThreeVector > > *  klmClusterInfo,
std::vector< std::map< const Track *, double > > *  bklmHitUsed 
)
private

Swim a single track (MUID) until it stops or leaves the target cylinder.

Definition at line 488 of file TrackExtrapolateG4e.cc.

492{
493 if (extState.pdgCode == 0)
494 return;
495 if (g4eState.GetMomentum().perp() <= m_MinPt)
496 return;
497 if (m_TargetMuid->GetDistanceFromPoint(g4eState.GetPosition()) < 0.0)
498 return;
499 G4ParticleDefinition* particle = G4ParticleTable::GetParticleTable()->FindParticle(extState.pdgCode);
500 double mass = particle->GetPDGMass();
501 double minPSq = (mass + m_MinKE) * (mass + m_MinKE) - mass * mass;
502 // Create structures for ECLCluster-Track matching
503 std::vector<double> eclClusterDistance;
504 std::vector<ExtHit> eclHit1, eclHit2, eclHit3;
505 if (eclClusterInfo != nullptr) {
506 eclClusterDistance.resize(eclClusterInfo->size(), 1.0E10); // "positive infinity"
507 ExtHit tempExtHit(extState.pdgCode, Const::EDetector::ECL, 0, EXT_FIRST,
508 extState.isCosmic, 0.0,
509 G4ThreeVector(), G4ThreeVector(), G4ErrorSymMatrix(6));
510 eclHit1.resize(eclClusterInfo->size(), tempExtHit);
511 eclHit2.resize(eclClusterInfo->size(), tempExtHit);
512 eclHit3.resize(eclClusterInfo->size(), tempExtHit);
513 }
514 // Create structures for KLMCluster-Track matching
515 std::vector<TrackClusterSeparation> klmHit;
516 if (klmClusterInfo != nullptr) {
517 klmHit.resize(klmClusterInfo->size()); // initialize each to huge distance
518 }
519 KLMMuidLikelihood* klmMuidLikelihood = m_klmMuidLikelihoods.appendNew(); // rest of this object will be filled later
520 klmMuidLikelihood->setPDGCode(extState.pdgCode);
521 if (extState.track != nullptr)
522 extState.track->addRelationTo(klmMuidLikelihood);
523 G4ErrorMode propagationMode = (extState.isCosmic ? G4ErrorMode_PropBackwards : G4ErrorMode_PropForwards);
524 m_ExtMgr->InitTrackPropagation(propagationMode);
525 while (true) {
526 const G4int errCode = m_ExtMgr->PropagateOneStep(&g4eState, propagationMode);
527 G4Track* track = g4eState.GetG4Track();
528 const G4Step* step = track->GetStep();
529 const G4StepPoint* preStepPoint = step->GetPreStepPoint();
530 const G4StepPoint* postStepPoint = step->GetPostStepPoint();
531 G4TouchableHandle preTouch = preStepPoint->GetTouchableHandle();
532 G4VPhysicalVolume* pVol = preTouch->GetVolume();
533 const G4int preStatus = preStepPoint->GetStepStatus();
534 const G4int postStatus = postStepPoint->GetStepStatus();
535 G4ThreeVector pos = track->GetPosition(); // this is at postStepPoint
536 G4ThreeVector mom = track->GetMomentum(); // ditto
537 // Ignore the zero-length step by PropagateOneStep() at each boundary
538 if (extState.isCosmic) mom = -mom;
539 if (step->GetStepLength() > 0.0) {
540 double dt = step->GetDeltaTime();
541 double dl = step->GetStepLength() / track->GetMaterial()->GetRadlen();
542 if (preStatus == fGeomBoundary) { // first step in this volume?
543 if (m_TargetExt->GetDistanceFromPoint(pos) < 0.0) { // only hits outside ECL during muid
544 if (m_EnterExit->find(pVol) != m_EnterExit->end()) {
545 createExtHit(EXT_ENTER, extState, g4eState, preStepPoint, preTouch);
546 }
547 }
548 }
549 if (extState.isCosmic) {
550 extState.tof -= dt;
551 extState.length -= dl;
552 } else {
553 extState.tof += dt;
554 extState.length += dl;
555 }
556 if (postStatus == fGeomBoundary) { // last step in this volume?
557 if (m_TargetExt->GetDistanceFromPoint(pos) < 0.0) { // only hits outside ECL during muid
558 if (m_EnterExit->find(pVol) != m_EnterExit->end()) {
559 createExtHit(EXT_EXIT, extState, g4eState, postStepPoint, preTouch);
560 }
561 }
562 }
563 if (createMuidHit(extState, g4eState, klmMuidLikelihood, bklmHitUsed)) {
564 // Force geant4e to update its G4Track from the Kalman-updated state
565 m_ExtMgr->GetPropagator()->SetStepN(0);
566 }
567 if (eclClusterInfo != nullptr) {
568 for (unsigned int c = 0; c < eclClusterInfo->size(); ++c) {
569 G4ThreeVector eclPos((*eclClusterInfo)[c].second);
570 G4ThreeVector prePos(preStepPoint->GetPosition());
571 G4ThreeVector diff(prePos - eclPos);
572 double distance = diff.mag();
573 if (distance < m_MaxECLTrackClusterDistance) {
574 // fallback ECLNEAR in case no ECLCROSS is found
575 if (distance < eclClusterDistance[c]) {
576 eclClusterDistance[c] = distance;
577 G4ErrorSymMatrix covariance(6, 0);
578 fromG4eToPhasespace(g4eState, covariance);
579 eclHit3[c].update(EXT_ECLNEAR, extState.tof, pos / CLHEP::cm, mom / CLHEP::GeV, covariance);
580 }
581 // find position of crossing of the track with the ECLCluster's sphere
582 if (eclHit1[c].getStatus() == EXT_FIRST) {
583 if (pos.mag2() >= eclPos.mag2()) {
584 double r = eclPos.mag();
585 double preD = prePos.mag() - r;
586 double postD = pos.mag() - r;
587 double f = postD / (postD - preD);
588 G4ThreeVector midPos = pos + (prePos - pos) * f;
589 double tof = extState.tof + dt * f * (extState.isCosmic ? +1 : -1); // in ns, at end of step
590 G4ErrorSymMatrix covariance(6, 0);
591 fromG4eToPhasespace(g4eState, covariance);
592 eclHit1[c].update(EXT_ECLCROSS, tof, midPos / CLHEP::cm, mom / CLHEP::GeV, covariance);
593 }
594 }
595 }
596 // find closest distance to the radial line to the ECLCluster
597 if (eclHit2[c].getStatus() == EXT_FIRST) {
598 G4ThreeVector delta(pos - prePos);
599 G4ThreeVector perp(eclPos.cross(delta));
600 double perpMag2 = perp.mag2();
601 if (perpMag2 > 1.0E-10) {
602 double dist = std::fabs(diff * perp) / std::sqrt(perpMag2);
603 if (dist < m_MaxECLTrackClusterDistance) {
604 double f = eclPos * (prePos.cross(perp)) / perpMag2;
605 if ((f > -0.5) && (f <= 1.0)) {
606 G4ThreeVector midPos(prePos + f * delta);
607 double length = extState.length + dl * (1.0 - f) * (extState.isCosmic ? +1 : -1);
608 G4ErrorSymMatrix covariance(6, 0);
609 fromG4eToPhasespace(g4eState, covariance);
610 eclHit2[c].update(EXT_ECLDL, length, midPos / CLHEP::cm, mom / CLHEP::GeV, covariance);
611 }
612 }
613 }
614 }
615 }
616 }
617 if (klmClusterInfo != nullptr) {
618 for (unsigned int c = 0; c < klmClusterInfo->size(); ++c) {
619 G4ThreeVector klmPos = (*klmClusterInfo)[c].second;
620 G4ThreeVector separation = klmPos - pos;
621 double distance = separation.mag();
622 if (distance < klmHit[c].getDistance()) {
623 klmHit[c].setDistance(distance);
624 klmHit[c].setTrackClusterAngle(mom.angle(separation));
625 klmHit[c].setTrackClusterSeparationAngle(mom.angle(klmPos));
626 klmHit[c].setTrackRotationAngle(extState.directionAtIP.angle(mom));
627 klmHit[c].setTrackClusterInitialSeparationAngle(extState.directionAtIP.angle(klmPos));
628 }
629 }
630 }
631 }
632 // Post-step momentum too low?
633 if (errCode || (mom.mag2() < minPSq)) {
634 break;
635 }
636 // Detect escapes from the imaginary target cylinder.
637 if (m_TargetMuid->GetDistanceFromPoint(pos) < 0.0) {
638 break;
639 }
640 // Stop extrapolating as soon as the track curls inward too much
641 if (pos.perp2() < m_MinRadiusSq) {
642 break;
643 }
644 } // track-extrapolation "infinite" loop
645
646 m_ExtMgr->EventTermination(propagationMode);
647
648 finishTrack(extState, klmMuidLikelihood, (g4eState.GetPosition().z() > m_OffsetZ));
649
650 if (eclClusterInfo != nullptr) {
651 for (unsigned int c = 0; c < eclClusterInfo->size(); ++c) {
652 if (eclHit1[c].getStatus() != EXT_FIRST) {
653 ExtHit* h = m_extHits.appendNew(eclHit1[c]);
654 (*eclClusterInfo)[c].first->addRelationTo(h);
655 if (extState.track != nullptr) {
656 extState.track->addRelationTo(h);
657 }
658 }
659 if (eclHit2[c].getStatus() != EXT_FIRST) {
660 ExtHit* h = m_extHits.appendNew(eclHit2[c]);
661 (*eclClusterInfo)[c].first->addRelationTo(h);
662 if (extState.track != nullptr) {
663 extState.track->addRelationTo(h);
664 }
665 }
666 if (eclHit3[c].getStatus() != EXT_FIRST) {
667 ExtHit* h = m_extHits.appendNew(eclHit3[c]);
668 (*eclClusterInfo)[c].first->addRelationTo(h);
669 if (extState.track != nullptr) {
670 extState.track->addRelationTo(h);
671 }
672 }
673 }
674 }
675
676 if (klmClusterInfo != nullptr) {
677 // here we set a relation only to the closest KLMCluster
678 // and we don't set any relation if the distance is too large
679 double minDistance = m_MaxKLMTrackClusterDistance;
680 unsigned int closestCluster = 0;
681 for (unsigned int c = 0; c < klmClusterInfo->size(); ++c) {
682 if (klmHit[c].getDistance() > 1.0E9) {
683 continue;
684 }
685 TrackClusterSeparation* h = m_trackClusterSeparations.appendNew(klmHit[c]);
686 (*klmClusterInfo)[c].first->addRelationTo(h); // relation KLMCluster to TrackSep
687 if (extState.track != nullptr) {
688 extState.track->addRelationTo(h); // relation Track to TrackSep
689 }
690 if (klmHit[c].getDistance() < minDistance) {
691 closestCluster = c;
692 minDistance = klmHit[c].getDistance();
693 }
694 }
695 if (minDistance < m_MaxKLMTrackClusterDistance) {
696 // set the relation Track to KLMCluster, using the distance as weight
697 // but for being consistent with the basf2 conventions, store the distance in cm
698 if (extState.track != nullptr) {
699 extState.track->addRelationTo((*klmClusterInfo)[closestCluster].first, 1. / (minDistance / CLHEP::cm));
700 }
701 }
702 }
703
704}
Class to store the likelihoods from KLM with additional information related to the extrapolation.
void setPDGCode(int pdg)
Set the PDG code of the particle hypothesis used during the extrapolation.
G4ErrorPropagator * GetPropagator() const
Get the propagator.
Definition: ExtManager.h:76
Store one Track-KLMCluster separation as a ROOT object.
bool createMuidHit(ExtState &, G4ErrorFreeTrajState &, KLMMuidLikelihood *, std::vector< std::map< const Track *, double > > *)
Create another MUID extrapolation hit for a track candidate.
void finishTrack(const ExtState &, KLMMuidLikelihood *, bool)
Complete muon identification after end of track extrapolation.

◆ terminate()

void terminate ( bool  flag)

Terminates this singleton.

Parameters
flagTrue if called by Muid module, false if called by Ext module.

Definition at line 417 of file TrackExtrapolateG4e.cc.

418{
419 if (m_DefaultHypotheses != nullptr)
420 delete m_DefaultHypotheses;
421 if (byMuid) {
422 delete m_TargetMuid;
423 for (auto const& muidBuilder : m_MuidBuilderMap)
424 delete muidBuilder.second;
425 }
426 if (m_TargetExt != nullptr) {
427 delete m_TargetExt;
428 m_TargetExt = nullptr;
429 }
430 if (m_EnterExit != nullptr) {
431 delete m_EnterExit;
432 delete m_BKLMVolumes;
434 m_EnterExit = nullptr;
435 m_BKLMVolumes = nullptr;
436 }
437}
void RunTermination()
Terminate a run and set state to G4ErrorState_Init.
Definition: ExtManager.cc:150

Member Data Documentation

◆ m_addHitsToRecoTrack

bool m_addHitsToRecoTrack = false
private

Parameter to add the found hits also to the reco tracks or not. Is turned off by default.

Definition at line 421 of file TrackExtrapolateG4e.h.

◆ m_BarrelHalfLength

double m_BarrelHalfLength
private

half-length (cm) of the barrel

Definition at line 372 of file TrackExtrapolateG4e.h.

◆ m_BarrelMaxR

double m_BarrelMaxR
private

maximum radius (cm) of the barrel

Definition at line 366 of file TrackExtrapolateG4e.h.

◆ m_BarrelMinR

double m_BarrelMinR
private

minimum radius (cm) of the barrel

Definition at line 369 of file TrackExtrapolateG4e.h.

◆ m_BarrelModuleMiddleRadius

double m_BarrelModuleMiddleRadius[2][BKLMElementNumbers::getMaximalSectorNumber()+1][BKLMElementNumbers::getMaximalLayerNumber()+1]
private

hit-plane radius (cm) at closest distance to IP of each barrel end | sector | layer

Definition at line 387 of file TrackExtrapolateG4e.h.

◆ m_BarrelNSector

int m_BarrelNSector
private

Number of barrel sectors.

Definition at line 363 of file TrackExtrapolateG4e.h.

◆ m_BarrelPhiStripVariance

double m_BarrelPhiStripVariance[BKLMElementNumbers::getMaximalLayerNumber()+1]
private

BKLM RPC phi-measuring strip position variance (cm^2) by layer.

Definition at line 378 of file TrackExtrapolateG4e.h.

◆ m_BarrelScintVariance

double m_BarrelScintVariance
private

BKLM scintillator strip position variance (cm^2)

Definition at line 384 of file TrackExtrapolateG4e.h.

◆ m_BarrelSectorPerp

G4ThreeVector m_BarrelSectorPerp[BKLMElementNumbers::getMaximalSectorNumber()+1]
private

normal unit vector of each barrel sector

Definition at line 391 of file TrackExtrapolateG4e.h.

◆ m_BarrelSectorPhi

G4ThreeVector m_BarrelSectorPhi[BKLMElementNumbers::getMaximalSectorNumber()+1]
private

azimuthal unit vector of each barrel sector

Definition at line 394 of file TrackExtrapolateG4e.h.

◆ m_BarrelZStripVariance

double m_BarrelZStripVariance[BKLMElementNumbers::getMaximalLayerNumber()+1]
private

BKLM RPC z-measuring strip position variance (cm^2) by layer.

Definition at line 381 of file TrackExtrapolateG4e.h.

◆ m_BeamPipeGeo

DBObjPtr<BeamPipeGeo> m_BeamPipeGeo
private

Conditions-database object for beam pipe geometry.

Definition at line 354 of file TrackExtrapolateG4e.h.

◆ m_BKLMVolumes

std::vector<G4VPhysicalVolume*>* m_BKLMVolumes
private

Pointers to BKLM geant4 sensitive (physical) volumes.

Definition at line 342 of file TrackExtrapolateG4e.h.

◆ m_COILGeometryPar

DBObjPtr<COILGeometryPar> m_COILGeometryPar
private

Conditions-database object for COIL geometry.

Definition at line 351 of file TrackExtrapolateG4e.h.

◆ m_DefaultHypotheses

std::vector<Const::ChargedStable>* m_DefaultHypotheses
private

Default ChargedStable hypotheses (needed as call argument but not used)

Definition at line 336 of file TrackExtrapolateG4e.h.

◆ m_eclClusters

StoreArray<ECLCluster> m_eclClusters
private

ECL clusters.

Definition at line 445 of file TrackExtrapolateG4e.h.

◆ m_eklmElementNumbers

const EKLMElementNumbers* m_eklmElementNumbers
private

EKLM element numbers.

Definition at line 427 of file TrackExtrapolateG4e.h.

◆ m_eklmTransformData

const EKLM::TransformDataGlobalAligned* m_eklmTransformData
private

EKLM transformation data.

Definition at line 433 of file TrackExtrapolateG4e.h.

◆ m_EndcapHalfLength

double m_EndcapHalfLength
private

half-length (cm) of either endcap

Definition at line 406 of file TrackExtrapolateG4e.h.

◆ m_EndcapMaxR

double m_EndcapMaxR
private

maximum radius (cm) of the endcaps

Definition at line 397 of file TrackExtrapolateG4e.h.

◆ m_EndcapMiddleZ

double m_EndcapMiddleZ
private

midpoint along z (cm) of the forward endcap from the KLM midpoint

Definition at line 403 of file TrackExtrapolateG4e.h.

◆ m_EndcapMinR

double m_EndcapMinR
private

minimum radius (cm) of the endcaps

Definition at line 400 of file TrackExtrapolateG4e.h.

◆ m_EndcapModuleMiddleZ

double m_EndcapModuleMiddleZ[BKLMElementNumbers::getMaximalLayerNumber()+1]
private

hit-plane z (cm) of each IP layer relative to KLM midpoint

Definition at line 418 of file TrackExtrapolateG4e.h.

◆ m_EndcapScintVariance

double m_EndcapScintVariance
private

EKLM scintillator strip position variance (cm^2)

Definition at line 415 of file TrackExtrapolateG4e.h.

◆ m_EnterExit

std::map<G4VPhysicalVolume*, enum VolTypes>* m_EnterExit
private

Pointers to geant4 physical volumes whose entry/exit points will be saved.

Definition at line 339 of file TrackExtrapolateG4e.h.

◆ m_extHits

StoreArray<ExtHit> m_extHits
private

Ext hits.

Definition at line 448 of file TrackExtrapolateG4e.h.

◆ m_ExtInitialized

bool m_ExtInitialized
private

Flag to indicate that EXT initialize() has been called.

Definition at line 297 of file TrackExtrapolateG4e.h.

◆ m_ExtMgr

Simulation::ExtManager* m_ExtMgr
private

Pointer to the ExtManager singleton.

Definition at line 327 of file TrackExtrapolateG4e.h.

◆ m_HypothesesExt

const std::vector<Const::ChargedStable>* m_HypothesesExt
private

ChargedStable hypotheses for EXT.

Definition at line 330 of file TrackExtrapolateG4e.h.

◆ m_HypothesesMuid

const std::vector<Const::ChargedStable>* m_HypothesesMuid
private

ChargedStable hypotheses for MUID.

Definition at line 333 of file TrackExtrapolateG4e.h.

◆ m_klmChannelStatus

DBObjPtr<KLMChannelStatus> m_klmChannelStatus
private

Conditions-database object for KLM channel status.

Definition at line 436 of file TrackExtrapolateG4e.h.

◆ m_klmClusters

StoreArray<KLMCluster> m_klmClusters
private

KLM clusters.

Definition at line 454 of file TrackExtrapolateG4e.h.

◆ m_klmElementNumbers

const KLMElementNumbers* m_klmElementNumbers
private

KLM element numbers.

Definition at line 430 of file TrackExtrapolateG4e.h.

◆ m_klmHit2ds

StoreArray<KLMHit2d> m_klmHit2ds
private

KLM 2d hits.

Definition at line 451 of file TrackExtrapolateG4e.h.

◆ m_klmLikelihoodParameters

DBObjPtr<KLMLikelihoodParameters> m_klmLikelihoodParameters
private

Conditions-database object for KLM likelihood parameters.

Definition at line 442 of file TrackExtrapolateG4e.h.

◆ m_klmMuidHits

StoreArray<KLMMuidHit> m_klmMuidHits
private

KLM muid hits.

Definition at line 457 of file TrackExtrapolateG4e.h.

◆ m_klmMuidLikelihoods

StoreArray<KLMMuidLikelihood> m_klmMuidLikelihoods
private

KLM muid likelihoods.

Definition at line 460 of file TrackExtrapolateG4e.h.

◆ m_klmStripEfficiency

DBObjPtr<KLMStripEfficiency> m_klmStripEfficiency
private

Conditions-database object for KLM strip efficiency.

Definition at line 439 of file TrackExtrapolateG4e.h.

◆ m_MagneticField

double m_MagneticField
private

Magnetic field z component (gauss) at origin.

Definition at line 309 of file TrackExtrapolateG4e.h.

◆ m_MaxDistSqInVariances

double m_MaxDistSqInVariances
private

user-defined maximum squared-distance (in number of variances) for matching hit to extrapolation

Definition at line 312 of file TrackExtrapolateG4e.h.

◆ m_MaxDt

double m_MaxDt
private

Coincidence window half-width for in-time KLM hits (ns)

Definition at line 306 of file TrackExtrapolateG4e.h.

◆ m_MaxECLTrackClusterDistance

double m_MaxECLTrackClusterDistance
private

user-defined maximum distance (mm) between ECLCluster and associated track (for EID)

Definition at line 318 of file TrackExtrapolateG4e.h.

◆ m_MaxKLMTrackClusterDistance

double m_MaxKLMTrackClusterDistance
private

user-defined maximum distance (mm) between KLMCluster and associated track (for KLID)

Definition at line 315 of file TrackExtrapolateG4e.h.

◆ m_MeanDt

double m_MeanDt
private

Mean hit - trigger time (ns)

Definition at line 303 of file TrackExtrapolateG4e.h.

◆ m_MinKE

double m_MinKE
private

Minimum kinetic energy in MeV for extrapolation to continue.

Definition at line 324 of file TrackExtrapolateG4e.h.

◆ m_MinPt

double m_MinPt
private

Minimum transverse momentum in MeV/c for extrapolation to be started.

Definition at line 321 of file TrackExtrapolateG4e.h.

◆ m_MinRadiusSq

double m_MinRadiusSq
private

Minimum squared radius (cm) outside of which extrapolation will continue.

Definition at line 357 of file TrackExtrapolateG4e.h.

◆ m_MuidBuilderMap

std::map<int, MuidBuilder*> m_MuidBuilderMap
private

PDF for the charged final state particle hypotheses.

Definition at line 424 of file TrackExtrapolateG4e.h.

◆ m_MuidInitialized

bool m_MuidInitialized
private

Flag to indicate that MUID initialize() has been called.

Definition at line 300 of file TrackExtrapolateG4e.h.

◆ m_OffsetZ

double m_OffsetZ
private

offset (cm) along z axis of KLM midpoint from IP

Definition at line 360 of file TrackExtrapolateG4e.h.

◆ m_OutermostActiveBackwardEndcapLayer

int m_OutermostActiveBackwardEndcapLayer
private

outermost backward-endcap layer that is active for muon identification (user-defined)

Definition at line 412 of file TrackExtrapolateG4e.h.

◆ m_OutermostActiveBarrelLayer

int m_OutermostActiveBarrelLayer
private

outermost barrel layer that is active for muon identification (user-defined)

Definition at line 375 of file TrackExtrapolateG4e.h.

◆ m_OutermostActiveForwardEndcapLayer

int m_OutermostActiveForwardEndcapLayer
private

outermost forward-endcap layer that is active for muon identification (user-defined)

Definition at line 409 of file TrackExtrapolateG4e.h.

◆ m_recoTracks

StoreArray<RecoTrack> m_recoTracks
private

Reco tracks.

Definition at line 463 of file TrackExtrapolateG4e.h.

◆ m_Singleton

TrackExtrapolateG4e * m_Singleton = nullptr
staticprivate

Stores pointer to the singleton class.

Definition at line 294 of file TrackExtrapolateG4e.h.

◆ m_TargetExt

Simulation::ExtCylSurfaceTarget* m_TargetExt
private

virtual "target" cylinder for EXT (boundary beyond which extrapolation ends)

Definition at line 345 of file TrackExtrapolateG4e.h.

◆ m_TargetMuid

Simulation::ExtCylSurfaceTarget* m_TargetMuid
private

virtual "target" cylinder for MUID (boundary beyond which extrapolation ends)

Definition at line 348 of file TrackExtrapolateG4e.h.

◆ m_trackClusterSeparations

StoreArray<TrackClusterSeparation> m_trackClusterSeparations
private

Track cluster separation.

Definition at line 469 of file TrackExtrapolateG4e.h.

◆ m_tracks

StoreArray<Track> m_tracks
private

Tracks.

Definition at line 466 of file TrackExtrapolateG4e.h.


The documentation for this class was generated from the following files: