File indexing completed on 2024-09-07 04:37:41
0001 #ifndef ME0Segment_ME0SegAlgoRU_h
0002 #define ME0Segment_ME0SegAlgoRU_h
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024 #include <RecoLocalMuon/GEMSegment/plugins/ME0SegmentAlgorithmBase.h>
0025 #include <DataFormats/GEMRecHit/interface/ME0RecHit.h>
0026 #include "MuonSegFit.h"
0027
0028 #include <vector>
0029
0030 class MuonSegFit;
0031 class ME0SegAlgoRU : public ME0SegmentAlgorithmBase {
0032 public:
0033 struct SegmentParameters {
0034 float maxETASeeds;
0035 float maxPhiSeeds;
0036 float maxPhiAdditional;
0037 float maxChi2Additional;
0038 float maxChi2Prune;
0039 float maxChi2GoodSeg;
0040 float maxTOFDiff;
0041 bool requireCentralBX;
0042 unsigned int minNumberOfHits;
0043 bool requireBeamConstr;
0044 };
0045
0046
0047 typedef std::vector<bool> BoolContainer;
0048 typedef std::vector<std::pair<float, HitAndPositionPtrContainer> > SegmentByMetricContainer;
0049
0050
0051 explicit ME0SegAlgoRU(const edm::ParameterSet& ps);
0052
0053 ~ME0SegAlgoRU() override {}
0054
0055
0056
0057
0058 std::vector<ME0Segment> run(const ME0Chamber* chamber, const HitAndPositionContainer& rechits) override;
0059
0060 private:
0061
0062 void lookForSegments(const SegmentParameters& params,
0063 const unsigned int n_seg_min,
0064 const HitAndPositionContainer& rechits,
0065 const std::vector<unsigned int>& recHits_per_layer,
0066 BoolContainer& used,
0067 std::vector<ME0Segment>& segments) const;
0068
0069 void tryAddingHitsToSegment(const float maxTOF,
0070 const float maxETA,
0071 const float maxPhi,
0072 const float maxChi2,
0073 std::unique_ptr<MuonSegFit>& current_fit,
0074 HitAndPositionPtrContainer& proto_segment,
0075 const BoolContainer& used,
0076 HitAndPositionContainer::const_iterator i1,
0077 HitAndPositionContainer::const_iterator i2) const;
0078
0079 void pruneBadHits(const float maxChi2,
0080 HitAndPositionPtrContainer& proto_segment,
0081 std::unique_ptr<MuonSegFit>& fit,
0082 const unsigned int n_seg_min) const;
0083
0084 void addUniqueSegments(SegmentByMetricContainer& proto_segments,
0085 std::vector<ME0Segment>& segments,
0086 BoolContainer& used) const;
0087
0088
0089 bool areHitsCloseInEta(const float maxETA, const bool beamConst, const GlobalPoint& h1, const GlobalPoint& h2) const;
0090 bool areHitsCloseInGlobalPhi(const float maxPHI,
0091 const unsigned int nLayDisp,
0092 const GlobalPoint& h1,
0093 const GlobalPoint& h2) const;
0094
0095
0096 std::unique_ptr<MuonSegFit> addHit(HitAndPositionPtrContainer& proto_segment, const HitAndPosition& aHit) const;
0097
0098 bool hasHitOnLayer(const HitAndPositionPtrContainer& proto_segment, const unsigned int layer) const;
0099
0100 std::unique_ptr<MuonSegFit> makeFit(const HitAndPositionPtrContainer& proto_segment) const;
0101
0102
0103 bool areHitsConsistentInTime(const float maxTOF,
0104 const HitAndPositionPtrContainer& proto_segment,
0105 const HitAndPosition& h) const;
0106
0107 bool isHitNearSegment(const float maxETA,
0108 const float maxPHI,
0109 const std::unique_ptr<MuonSegFit>& fit,
0110 const HitAndPositionPtrContainer& proto_segment,
0111 const HitAndPosition& h) const;
0112
0113 float getHitSegChi2(const std::unique_ptr<MuonSegFit>& fit, const ME0RecHit& hit) const;
0114
0115 GlobalPoint globalAtZ(const std::unique_ptr<MuonSegFit>& fit, float z) const;
0116
0117
0118 void compareProtoSegment(std::unique_ptr<MuonSegFit>& current_fit,
0119 HitAndPositionPtrContainer& current_proto_segment,
0120 const HitAndPosition& new_hit) const;
0121
0122 void increaseProtoSegment(const float maxChi2,
0123 std::unique_ptr<MuonSegFit>& current_fit,
0124 HitAndPositionPtrContainer& current_proto_segment,
0125 const HitAndPosition& new_hit) const;
0126
0127 const std::string myName;
0128 bool doCollisions;
0129 bool allowWideSegments;
0130
0131 SegmentParameters stdParameters;
0132 SegmentParameters displacedParameters;
0133 SegmentParameters wideParameters;
0134
0135
0136 const ME0Chamber* theChamber;
0137 };
0138
0139 #endif