ATLAS Offline Software
Loading...
Searching...
No Matches
FastCaloSimCaloExtrapolation.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2023 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef FastCaloSimCaloExtrapolation_H
6#define FastCaloSimCaloExtrapolation_H
7
9
14
15
17
18class TFCSTruthState;
19class G4FieldTrack;
20
21
29
30
31class FastCaloSimCaloExtrapolation: public extends<AthAlgTool, IFastCaloSimCaloExtrapolation>
32{
33
34public:
35
36 FastCaloSimCaloExtrapolation(const std::string& t, const std::string& n, const IInterface* p);
38
39 virtual StatusCode initialize() override final;
40 virtual StatusCode finalize() override final;
41
42 enum SUBPOS {
43 SUBPOS_MID = TFCSExtrapolationState::SUBPOS_MID, //MID=middle of calo layer
44 SUBPOS_ENT = TFCSExtrapolationState::SUBPOS_ENT, //ENT=entrance of calo layer
46 };
47
49 INSIDE, //hit position is inside cylinder bounds
50 OUTSIDE, //hit position is outside cylinder bounds
51 ON //hit position is on cylinder bounds
52 };
53
54 virtual void extrapolate(TFCSExtrapolationState& result, const TFCSTruthState* truth, const std::vector<G4FieldTrack>& caloSteps) const override final;
55 virtual void extrapolate(TFCSExtrapolationState& result, const TFCSTruthState* truth) const override final;
56
57protected:
58
59 const IFastCaloSimGeometryHelper* GetCaloGeometry() const {return &(*m_CaloGeometryHelper);};
60
61 /*Main extrapolation methods*/
62
64 bool extrapolateToCylinder(const std::vector<G4FieldTrack>& caloSteps, float cylR, float cylZ, Amg::Vector3D& extPos, Amg::Vector3D& momDir) const;
66 void extrapolateToID(TFCSExtrapolationState& result, const std::vector<G4FieldTrack>& caloSteps, const TFCSTruthState* truth) const;
68 void extrapolateToLayers(TFCSExtrapolationState& result, const std::vector<G4FieldTrack>& caloSteps, const TFCSTruthState* truth) const;
69
70 /*Extrapolator helper methods*/
71
73 void findPCA(float cylR, float cylZ, Amg::Vector3D& hitPos1, Amg::Vector3D& hitPos2, Amg::Vector3D& PCA) const;
75 static double getPointLineSegmentDistance(Amg::Vector3D& point, Amg::Vector3D& hitPos1, Amg::Vector3D& hitPos2) ;
77 void getIterativePCA(float cylR, float cylZ, Amg::Vector3D& BoundA, Amg::Vector3D& BoundB, Amg::Vector3D& PCA) const;
79 static bool isOnSegment(Amg::Vector3D& point, Amg::Vector3D& hitPos1, Amg::Vector3D& hitPos2) ;
81 static bool cylinderEndcapIntersection(float cylR, float cylZ, bool positiveEndcap, Amg::Vector3D& pointA, Amg::Vector3D& pointB, Amg::Vector3D& intersection) ;
84 bool extrapolateWithIntersection(const std::vector<G4FieldTrack>& caloSteps, float cylR, float cylZ, Amg::Vector3D& extPos, Amg::Vector3D& momDir) const;
86 bool extrapolateWithPCA(const std::vector<G4FieldTrack>& caloSteps, float cylR, float cylZ, Amg::Vector3D& extPos, Amg::Vector3D& momDir) const;
88 static bool doesTravelThroughSurface(float cylR, float cylZ, Amg::Vector3D& hitPos1, Amg::Vector3D& hitPos2) ;
90 int whichIntersection(float cylR, float cylZ, Amg::Vector3D& hitPos1, Amg::Vector3D& hitPos2, Amg::Vector3D& intersectionA, Amg::Vector3D intersectionB) const;
92 int circleLineIntersection2D(float circR, Amg::Vector3D& pointA, Amg::Vector3D& pointB, Amg::Vector3D& intersectA, Amg::Vector3D& intersectB) const;
94 int cylinderLineIntersection(float cylR, float cylZ, Amg::Vector3D& pointA, Amg::Vector3D& pointB, Amg::Vector3D& intersectA, Amg::Vector3D& intersectB) const;
96 static enum HITPOSITION whereOnCylinder(float cylR, float cylZ, Amg::Vector3D& hitPos) ;
98 static Amg::Vector3D projectOnCylinder(float cylR, float cylZ, Amg::Vector3D& hitPos) ;
100 CylinderIntersections getCylinderIntersections(float cylR, float cylZ, Amg::Vector3D& hitPos1, Amg::Vector3D& hitPos2) const;
101
102 //helper methods for calo geometry
103 void minmaxeta(int sample, double eta, double& mineta, double& maxeta) const;
104 bool isCaloBarrel(int sample) const;
105 double deta (int sample, double eta) const;
106 double rzmid(int sample, double eta) const;
107 double rzent(int sample, double eta) const;
108 double rzext(int sample, double eta) const;
109 double rmid (int sample, double eta) const;
110 double rent (int sample, double eta) const;
111 double rext (int sample, double eta) const;
112 double zmid (int sample, double eta) const;
113 double zent (int sample, double eta) const;
114 double zext (int sample, double eta) const;
115 double rpos (int sample, double eta, int subpos = CaloSubPos::SUBPOS_MID) const;
116 double zpos (int sample, double eta, int subpos = CaloSubPos::SUBPOS_MID) const;
117 double rzpos(int sample, double eta, int subpos = CaloSubPos::SUBPOS_MID) const;
118
119 HepPDT::ParticleDataTable* m_particleDataTable{nullptr};
120
121 //uniquely defined ID-Calo surfaces
122 FloatArrayProperty m_CaloBoundaryR{this, "CaloBoundaryR", {1148.0,120.0,41.0}};
123 FloatArrayProperty m_CaloBoundaryZ{this, "CaloBoundaryZ", {3550.0,4587.0,4587.0}};
124
125 // FastCaloSim particle transport with ATLAS tracking tools
126 PublicToolHandle<IFastCaloSimCaloTransportation> m_CaloTransportation{this, "CaloTransportation", "FastCaloSimCaloTransportation"};
127 // FastCaloSim geometry helper
128 PublicToolHandle<IFastCaloSimGeometryHelper> m_CaloGeometryHelper{this, "CaloGeometryHelper", "FastCaloSimGeometryHelper"};
129
130};
131
132#endif // FastCaloSimCaloExtrapolation_H
Scalar eta() const
pseudorapidity method
Cached pointer with atomic update.
virtual void extrapolate(TFCSExtrapolationState &result, const TFCSTruthState *truth, const std::vector< G4FieldTrack > &caloSteps) const override final
double zpos(int sample, double eta, int subpos=CaloSubPos::SUBPOS_MID) const
static enum HITPOSITION whereOnCylinder(float cylR, float cylZ, Amg::Vector3D &hitPos)
Checks if position of hitPos is inside, outside or on the cylinder bounds.
bool extrapolateWithPCA(const std::vector< G4FieldTrack > &caloSteps, float cylR, float cylZ, Amg::Vector3D &extPos, Amg::Vector3D &momDir) const
Extrapolates to the cylinder using the PCA to the polygon spanned by the individual line segments fro...
void extrapolateToID(TFCSExtrapolationState &result, const std::vector< G4FieldTrack > &caloSteps, const TFCSTruthState *truth) const
Extrapolates to ID using three uniquely defined cylinder surfaces.
double deta(int sample, double eta) const
double rext(int sample, double eta) const
void getIterativePCA(float cylR, float cylZ, Amg::Vector3D &BoundA, Amg::Vector3D &BoundB, Amg::Vector3D &PCA) const
Finds PCA iteratively given two bounds A and B on a line segment, used for (rare) cases with no easy ...
void extrapolateToLayers(TFCSExtrapolationState &result, const std::vector< G4FieldTrack > &caloSteps, const TFCSTruthState *truth) const
Extrapolates to all other layers of the calorimeter.
static Amg::Vector3D projectOnCylinder(float cylR, float cylZ, Amg::Vector3D &hitPos)
Projects position hitPos onto the cylinder surface and returns projected position.
virtual StatusCode finalize() override final
CylinderIntersections getCylinderIntersections(float cylR, float cylZ, Amg::Vector3D &hitPos1, Amg::Vector3D &hitPos2) const
Analytically computes the intersection between the (infinite) line spanned by hitPos1 and hitPos2 wit...
double rpos(int sample, double eta, int subpos=CaloSubPos::SUBPOS_MID) const
void minmaxeta(int sample, double eta, double &mineta, double &maxeta) const
HepPDT::ParticleDataTable * m_particleDataTable
const IFastCaloSimGeometryHelper * GetCaloGeometry() const
FastCaloSimCaloExtrapolation(const std::string &t, const std::string &n, const IInterface *p)
PublicToolHandle< IFastCaloSimGeometryHelper > m_CaloGeometryHelper
double zent(int sample, double eta) const
int circleLineIntersection2D(float circR, Amg::Vector3D &pointA, Amg::Vector3D &pointB, Amg::Vector3D &intersectA, Amg::Vector3D &intersectB) const
Analytically computes 2D intersections between circle of radius circR and (infinite) line spanned by ...
bool extrapolateWithIntersection(const std::vector< G4FieldTrack > &caloSteps, float cylR, float cylZ, Amg::Vector3D &extPos, Amg::Vector3D &momDir) const
Extrapolates position on cylinder by finding intersections of subsequent hit positions,...
double rzmid(int sample, double eta) const
double rzpos(int sample, double eta, int subpos=CaloSubPos::SUBPOS_MID) const
int whichIntersection(float cylR, float cylZ, Amg::Vector3D &hitPos1, Amg::Vector3D &hitPos2, Amg::Vector3D &intersectionA, Amg::Vector3D intersectionB) const
Returns ID of more sensible intersection between line segment spanned by hitPos1 and hitPos2 and cyli...
void findPCA(float cylR, float cylZ, Amg::Vector3D &hitPos1, Amg::Vector3D &hitPos2, Amg::Vector3D &PCA) const
Finds Point of Closest Approach (PCA) on the cylinder defined by radius cylR and half-length cylZ of ...
double rent(int sample, double eta) const
virtual StatusCode initialize() override final
bool extrapolateToCylinder(const std::vector< G4FieldTrack > &caloSteps, float cylR, float cylZ, Amg::Vector3D &extPos, Amg::Vector3D &momDir) const
Finds best extrapolation extPos from the caloSteps for a cylinder defined by radius cylR and half-len...
static double getPointLineSegmentDistance(Amg::Vector3D &point, Amg::Vector3D &hitPos1, Amg::Vector3D &hitPos2)
Computes the distance between a point and the line segment spanned by hitPos1 and hitPos2.
int cylinderLineIntersection(float cylR, float cylZ, Amg::Vector3D &pointA, Amg::Vector3D &pointB, Amg::Vector3D &intersectA, Amg::Vector3D &intersectB) const
Analytically computes the intersection between the (infinite) line defined by pointA and pointB and t...
static bool doesTravelThroughSurface(float cylR, float cylZ, Amg::Vector3D &hitPos1, Amg::Vector3D &hitPos2)
Returns true if the line segment spanned by hitPos1 and hitPos2 crosses the cylinder surface,...
double rzent(int sample, double eta) const
double zext(int sample, double eta) const
double rzext(int sample, double eta) const
double rmid(int sample, double eta) const
PublicToolHandle< IFastCaloSimCaloTransportation > m_CaloTransportation
double zmid(int sample, double eta) const
static bool cylinderEndcapIntersection(float cylR, float cylZ, bool positiveEndcap, Amg::Vector3D &pointA, Amg::Vector3D &pointB, Amg::Vector3D &intersection)
Computes intersection between the (infinite) line spanned by pointA and pointB with the positive (neg...
static bool isOnSegment(Amg::Vector3D &point, Amg::Vector3D &hitPos1, Amg::Vector3D &hitPos2)
Returns true if point lies on the line segment spanned by hitPos1 and hitPos2, otherwise returns fals...
std::vector< std::string > intersection(std::vector< std::string > &v1, std::vector< std::string > &v2)
Eigen::Matrix< double, 3, 1 > Vector3D