CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TotemTimingLocalTrackFitter.cc
Go to the documentation of this file.
1 /****************************************************************************
2  *
3  * This is a part of CTPPS offline software.
4  * Authors:
5  * Laurent Forthomme (laurent.forthomme@cern.ch)
6  * Nicola Minafra (nicola.minafra@cern.ch)
7  * Mateusz Szpyrka (mateusz.szpyrka@cern.ch)
8  * Christopher Misan (krzysztof.misan@cern.ch)
9  *
10  ****************************************************************************/
11 
12 #include <memory>
13 #include <utility>
14 
19 
22 
24 
29 
31 template <typename T>
33 public:
35 
37 
38 private:
40  void produce(edm::Event&, const edm::EventSetup&) override;
41 
45  std::map<CTPPSDetId, TotemTimingTrackRecognition> trk_algo_map_;
46 };
47 template <typename T>
49  : recHitsToken_(consumes<edm::DetSetVector<TotemTimingRecHit> >(iConfig.getParameter<edm::InputTag>("recHitsTag"))),
50  maxPlaneActiveChannels_(iConfig.getParameter<int>("maxPlaneActiveChannels")),
51  trk_algo_params_(iConfig.getParameter<edm::ParameterSet>("trackingAlgorithmParams")) {
52  produces<edm::DetSetVector<TotemTimingLocalTrack> >();
53 }
54 template <typename T>
56  auto pOut = std::make_unique<edm::DetSetVector<TotemTimingLocalTrack> >();
57 
59  iEvent.getByToken(recHitsToken_, recHits);
60 
61  std::map<T, int> planeActivityMap;
62 
63  auto motherId = [](const edm::det_id_type& detid) {
64  T out(detid);
65  out.setChannel(0);
66  return out;
67  };
68 
69  for (const auto& vec : *recHits)
70  planeActivityMap[motherId(vec.detId())] += vec.size();
71 
72  // feed hits to the track producers
73  for (const auto& vec : *recHits) {
74  const CTPPSDetId raw_detid(vec.detId());
75  T detid(raw_detid.arm(), raw_detid.station(), raw_detid.rp());
76  // if algorithm is not found, build it
77  if (trk_algo_map_.count(detid) == 0)
78  trk_algo_map_.insert(std::make_pair(detid, trk_algo_params_));
79 
80  auto detId = motherId(vec.detId());
81  if (planeActivityMap[detId] > maxPlaneActiveChannels_)
82  continue;
83 
84  detId.setPlane(0);
85  for (const auto& hit : vec) {
86  if (trk_algo_map_.find(detId) == trk_algo_map_.end())
87  throw cms::Exception("TotemTimingLocalTrackFitter")
88  << "Invalid detId for rechit: arm=" << detId.arm() << ", rp=" << detId.rp();
89  trk_algo_map_.find(detId)->second.addHit(hit);
90  }
91  }
92 
93  // retrieves tracks for all hit sets
94  for (auto& trk_algo_entry : trk_algo_map_) {
95  pOut->find_or_insert(trk_algo_entry.first);
96  trk_algo_entry.second.produceTracks(pOut->operator[](trk_algo_entry.first));
97  }
98 
99  iEvent.put(std::move(pOut));
100 
101  // remove all hits from the track producers to prepare for the next event
102  for (auto& trk_algo_entry : trk_algo_map_)
103  trk_algo_entry.second.clear();
104 }
105 template <typename T>
109  desc.add<edm::InputTag>("recHitsTag", edm::InputTag("totemTimingRecHits"))
110  ->setComment("input rechits collection to retrieve");
111  desc.add<int>("maxPlaneActiveChannels", 2)->setComment("threshold for discriminating noisy planes");
112 
113  edm::ParameterSetDescription trackingAlgoParams;
114  trackingAlgoParams.add<double>("threshold", 1.5)
115  ->setComment("minimal number of rechits to be observed before launching the track recognition algorithm");
116  trackingAlgoParams.add<double>("thresholdFromMaximum", 0.5)
117  ->setComment("threshold relative to hit profile function local maximum for determining the width of the track");
118  trackingAlgoParams.add<double>("resolution", 0.01 /* mm */)
119  ->setComment("spatial resolution on the horizontal coordinate (in mm)");
120  trackingAlgoParams.add<double>("sigma", 0.)
121  ->setComment("pixel efficiency function parameter determining the smoothness of the step");
122  trackingAlgoParams.add<double>("tolerance", 0.1 /* mm */)
123  ->setComment("tolerance used for checking if the track contains certain hit");
124 
125  trackingAlgoParams.add<std::string>("pixelEfficiencyFunction", "(x>[0]-0.5*[1]-0.05)*(x<[0]+0.5*[1]-0.05)+0*[2]")
126  ->setComment(
127  "efficiency function for single pixel\n"
128  "can be defined as:\n"
129  " * Precise: "
130  "(TMath::Erf((x-[0]+0.5*([1]-0.05))/([2]/4)+2)+1)*TMath::Erfc((x-[0]-0.5*([1]-0.05))/([2]/4)-2)/4\n"
131  " * Fast: "
132  "(x>[0]-0.5*([1]-0.05))*(x<[0]+0.5*([1]-0.05))+((x-[0]+0.5*([1]-0.05)+[2])/"
133  "[2])*(x>[0]-0.5*([1]-0.05)-[2])*(x<[0]-0.5*([1]-0.05))+(2-(x-[0]-0.5*([1]-0.05)+[2])/"
134  "[2])*(x>[0]+0.5*([1]-0.05))*(x<[0]+0.5*([1]-0.05)+[2])\n"
135  " * Legacy: (1/(1+exp(-(x-[0]+0.5*([1]-0.05))/[2])))*(1/(1+exp((x-[0]-0.5*([1]-0.05))/[2])))\n"
136  " * Default (sigma ignored): (x>[0]-0.5*[1]-0.05)*(x<[0]+0.5*[1]-0.05)+0*[2]\n"
137  "with:\n"
138  " [0]: centre of pad\n"
139  " [1]: width of pad\n"
140  " [2]: sigma: distance between efficiency ~100 -> 0 outside width");
141 
142  trackingAlgoParams.add<double>("yPosition", 0.0)->setComment("vertical offset of the outcoming track centre");
143  trackingAlgoParams.add<double>("yWidth", 0.0)->setComment("vertical track width");
144  desc.add<edm::ParameterSetDescription>("trackingAlgorithmParams", trackingAlgoParams)
145  ->setComment("list of parameters associated to the track recognition algorithm");
146  return desc;
147 }
148 template <>
150  auto desc = fillDescriptionsShared(descr);
151  descr.add("totemTimingLocalTracks", desc);
152 }
153 
154 template <>
156  auto desc = fillDescriptionsShared(descr);
157  descr.add("diamondSampicLocalTracks", desc);
158 }
void setComment(std::string const &value)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
static edm::ParameterSetDescription fillDescriptionsShared(edm::ConfigurationDescriptions &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
int iEvent
Definition: GenABIO.cc:224
def move
Definition: eostools.py:511
ParameterDescriptionBase * add(U const &iLabel, T const &value)
TotemTimingLocalTrackFitter(const edm::ParameterSet &)
uint32_t det_id_type
Definition: DetSet.h:20
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const edm::ParameterSet trk_algo_params_
Base class for CTPPS detector IDs.
Definition: CTPPSDetId.h:32
void produce(edm::Event &, const edm::EventSetup &) override
long double T
const edm::EDGetTokenT< edm::DetSetVector< TotemTimingRecHit > > recHitsToken_
static void fillDescriptions(edm::ConfigurationDescriptions &)
std::map< CTPPSDetId, TotemTimingTrackRecognition > trk_algo_map_