CMS 3D CMS Logo

CTPPSDiamondLocalTrackFitter.cc
Go to the documentation of this file.
1 /****************************************************************************
2  *
3  * This is a part of PPS offline software.
4  * Authors:
5  * Laurent Forthomme (laurent.forthomme@cern.ch)
6  * Nicola Minafra (nicola.minafra@cern.ch)
7  *
8  ****************************************************************************/
9 
10 #include <memory>
11 
16 
19 
21 
25 
27 
29 public:
31 
33 
34 private:
35  void produce(edm::Event&, const edm::EventSetup&) override;
36 
39  std::unordered_map<CTPPSDetId, std::unique_ptr<CTPPSDiamondTrackRecognition> > trk_algo_;
40 };
41 
43  : recHitsToken_(
44  consumes<edm::DetSetVector<CTPPSDiamondRecHit> >(iConfig.getParameter<edm::InputTag>("recHitsTag"))),
45  trk_algo_params_(iConfig.getParameter<edm::ParameterSet>("trackingAlgorithmParams")) {
46  produces<edm::DetSetVector<CTPPSDiamondLocalTrack> >();
47 }
48 
50  // prepare the output
51  auto pOut = std::make_unique<edm::DetSetVector<CTPPSDiamondLocalTrack> >();
52 
54  iEvent.getByToken(recHitsToken_, recHits);
55 
56  // clear all hits possibly inherited from previous event
57  for (auto& algo_vs_id : trk_algo_)
58  algo_vs_id.second->clear();
59 
60  // feed hits to the track producers
61  for (const auto& vec : *recHits) {
62  const CTPPSDiamondDetId raw_detid(vec.detId()), detid(raw_detid.arm(), raw_detid.station(), raw_detid.rp());
63  // if algorithm is not found, build it
64  if (trk_algo_.count(detid) == 0)
65  trk_algo_[detid] = std::make_unique<CTPPSDiamondTrackRecognition>(trk_algo_params_);
66  for (const auto& hit : vec)
67  // skip hits without a leading edge
69  trk_algo_[detid]->addHit(hit);
70  }
71 
72  // build the tracks for all stations
73  for (auto& algo_vs_id : trk_algo_) {
74  auto& tracks = pOut->find_or_insert(algo_vs_id.first);
75  algo_vs_id.second->produceTracks(tracks);
76  }
77 
78  iEvent.put(std::move(pOut));
79 }
80 
83  desc.add<edm::InputTag>("recHitsTag", edm::InputTag("ctppsDiamondRecHits"))
84  ->setComment("input rechits collection to retrieve");
85 
86  edm::ParameterSetDescription trackingAlgoParams;
87  trackingAlgoParams.add<double>("threshold", 1.5)
88  ->setComment("minimal number of rechits to be observed before launching the track recognition algorithm");
89  trackingAlgoParams.add<double>("thresholdFromMaximum", 0.5);
90  trackingAlgoParams.add<double>("resolution", 0.01 /* mm */)
91  ->setComment("spatial resolution on the horizontal coordinate (in mm)");
92  trackingAlgoParams.add<double>("sigma", 0.1);
93  trackingAlgoParams.add<double>("startFromX", -0.5 /* mm */)
94  ->setComment("starting horizontal coordinate of rechits for the track recognition");
95  trackingAlgoParams.add<double>("stopAtX", 19.5 /* mm */)
96  ->setComment("ending horizontal coordinate of rechits for the track recognition");
97  trackingAlgoParams.add<double>("tolerance", 0.1 /* mm */)
98  ->setComment("tolerance used for checking if the track contains certain hit");
99 
100  trackingAlgoParams.add<std::string>("pixelEfficiencyFunction", "(x>[0]-0.5*[1])*(x<[0]+0.5*[1])+0*[2]")
101  ->setComment(
102  "efficiency function for single pixel\n"
103  "can be defined as:\n"
104  " * Precise: (TMath::Erf((x-[0]+0.5*[1])/([2]/4)+2)+1)*TMath::Erfc((x-[0]-0.5*[1])/([2]/4)-2)/4\n"
105  " * Fast: "
106  "(x>[0]-0.5*[1])*(x<[0]+0.5*[1])+((x-[0]+0.5*[1]+[2])/"
107  "[2])*(x>[0]-0.5*[1]-[2])*(x<[0]-0.5*[1])+(2-(x-[0]-0.5*[1]+[2])/[2])*(x>[0]+0.5*[1])*(x<[0]+0.5*[1]+[2])\n"
108  " * Legacy: (1/(1+exp(-(x-[0]+0.5*[1])/[2])))*(1/(1+exp((x-[0]-0.5*[1])/[2])))\n"
109  " * Default (sigma ignored): (x>[0]-0.5*[1])*(x<[0]+0.5*[1])+0*[2]\n"
110  "with:\n"
111  " [0]: centre of pad\n"
112  " [1]: width of pad\n"
113  " [2]: sigma: distance between efficiency ~100 -> 0 outside width");
114 
115  trackingAlgoParams.add<double>("yPosition", 0.0)->setComment("vertical offset of the outcoming track centre");
116  trackingAlgoParams.add<double>("yWidth", 0.0)->setComment("vertical track width");
117  trackingAlgoParams.add<bool>("excludeSingleEdgeHits", true)
118  ->setComment("exclude rechits with missing leading/trailing edge");
119 
120  desc.add<edm::ParameterSetDescription>("trackingAlgorithmParams", trackingAlgoParams)
121  ->setComment("list of parameters associated to the track recognition algorithm");
122 
123  descr.add("ctppsDiamondLocalTracks", desc);
124 }
125 
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
CTPPSDiamondTrackRecognition.h
edm::EDGetTokenT
Definition: EDGetToken.h:33
CTPPSDiamondRecHit
Reconstructed hit in diamond detectors.
Definition: CTPPSDiamondRecHit.h:16
edm
HLT enums.
Definition: AlignableModifier.h:19
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89281
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
EDProducer.h
CTPPSDiamondLocalTrackFitter::CTPPSDiamondLocalTrackFitter
CTPPSDiamondLocalTrackFitter(const edm::ParameterSet &)
Definition: CTPPSDiamondLocalTrackFitter.cc:42
edm::Handle
Definition: AssociativeIterator.h:50
CTPPSDiamondLocalTrackFitter::trk_algo_
std::unordered_map< CTPPSDetId, std::unique_ptr< CTPPSDiamondTrackRecognition > > trk_algo_
Definition: CTPPSDiamondLocalTrackFitter.cc:45
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::ConfigurationDescriptions::add
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: ConfigurationDescriptions.cc:57
CTPPSDiamondRecHit.h
edm::ConfigurationDescriptions
Definition: ConfigurationDescriptions.h:28
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
FastTrackerRecHitMaskProducer_cfi.recHits
recHits
Definition: FastTrackerRecHitMaskProducer_cfi.py:8
CTPPSDiamondLocalTrackFitter
Definition: CTPPSDiamondLocalTrackFitter.cc:28
CTPPSDiamondDetId.h
CTPPSDiamondDetId
Detector ID class for CTPPS Timing Diamond detectors. Bits [19:31] : Assigend in CTPPSDetId Calss Bit...
Definition: CTPPSDiamondDetId.h:24
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
tracks
const uint32_t *__restrict__ const HitContainer *__restrict__ TkSoA *__restrict__ tracks
Definition: CAHitNtupletGeneratorKernelsImpl.h:159
ParameterSet
Definition: Functions.h:16
CTPPSDiamondRecHit::TIMESLICE_WITHOUT_LEADING
static constexpr int TIMESLICE_WITHOUT_LEADING
Definition: CTPPSDiamondRecHit.h:44
CTPPSDiamondLocalTrack.h
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::stream::EDProducer
Definition: EDProducer.h:38
edm::EventSetup
Definition: EventSetup.h:58
DetSetVector.h
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
CTPPSDiamondLocalTrackFitter::trk_algo_params_
const edm::ParameterSet trk_algo_params_
Definition: CTPPSDiamondLocalTrackFitter.cc:44
CTPPSDiamondLocalTrackFitter::fillDescriptions
static void fillDescriptions(edm::ConfigurationDescriptions &)
Definition: CTPPSDiamondLocalTrackFitter.cc:81
CTPPSDiamondLocalTrackFitter::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: CTPPSDiamondLocalTrackFitter.cc:49
Frameworkfwd.h
CTPPSDiamondLocalTrackFitter::recHitsToken_
edm::EDGetTokenT< edm::DetSetVector< CTPPSDiamondRecHit > > recHitsToken_
Definition: CTPPSDiamondLocalTrackFitter.cc:43
ParameterSet.h
edm::ParameterDescriptionNode::setComment
void setComment(std::string const &value)
Definition: ParameterDescriptionNode.cc:106
edm::Event
Definition: Event.h:73
StreamID.h
edm::InputTag
Definition: InputTag.h:15
hit
Definition: SiStripHitEffFromCalibTree.cc:88