CMS 3D CMS Logo

MkFitProducer.cc
Go to the documentation of this file.
2 
7 
12 
19 
20 // mkFit includes
25 
26 // TBB includes
27 #include "oneapi/tbb/task_arena.h"
28 
29 // std includes
30 #include <functional>
31 
32 class MkFitProducer : public edm::global::EDProducer<edm::StreamCache<mkfit::MkBuilderWrapper>> {
33 public:
34  explicit MkFitProducer(edm::ParameterSet const& iConfig);
35  ~MkFitProducer() override;
36 
37  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
38 
39  std::unique_ptr<mkfit::MkBuilderWrapper> beginStream(edm::StreamID) const override;
40 
41 private:
42  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
43 
44  void stripClusterChargeCut(const std::vector<float>& stripClusterCharge, std::vector<bool>& mask) const;
45 
56  const float minGoodStripCharge_;
57  const bool seedCleaning_;
58  const bool backwardFitInCMSSW_;
59  const bool removeDuplicates_;
60  const bool mkFitSilent_;
61  const bool limitConcurrency_;
62 };
63 
65  : pixelHitsToken_{consumes(iConfig.getParameter<edm::InputTag>("pixelHits"))},
66  stripHitsToken_{consumes(iConfig.getParameter<edm::InputTag>("stripHits"))},
67  stripClusterChargeToken_{consumes(iConfig.getParameter<edm::InputTag>("stripHits"))},
68  eventOfHitsToken_{consumes(iConfig.getParameter<edm::InputTag>("eventOfHits"))},
69  seedToken_{consumes(iConfig.getParameter<edm::InputTag>("seeds"))},
70  mkFitGeomToken_{esConsumes()},
71  mkFitIterConfigToken_{esConsumes(iConfig.getParameter<edm::ESInputTag>("config"))},
72  putToken_{produces<MkFitOutputWrapper>()},
73  minGoodStripCharge_{static_cast<float>(
74  iConfig.getParameter<edm::ParameterSet>("minGoodStripCharge").getParameter<double>("value"))},
75  seedCleaning_{iConfig.getParameter<bool>("seedCleaning")},
76  backwardFitInCMSSW_{iConfig.getParameter<bool>("backwardFitInCMSSW")},
77  removeDuplicates_{iConfig.getParameter<bool>("removeDuplicates")},
78  mkFitSilent_{iConfig.getUntrackedParameter<bool>("mkFitSilent")},
79  limitConcurrency_{iConfig.getUntrackedParameter<bool>("limitConcurrency")} {
80  const auto clustersToSkip = iConfig.getParameter<edm::InputTag>("clustersToSkip");
81  if (not clustersToSkip.label().empty()) {
82  pixelMaskToken_ = consumes(clustersToSkip);
83  stripMaskToken_ = consumes(clustersToSkip);
84  }
85 
86  const auto build = iConfig.getParameter<std::string>("buildingRoutine");
87  if (build == "bestHit") {
88  //buildFunction_ = mkfit::runBuildingTestPlexBestHit;
89  throw cms::Exception("Configuration") << "bestHit is temporarily disabled";
90  } else if (build == "standard") {
91  //buildFunction_ = mkfit::runBuildingTestPlexStandard;
92  throw cms::Exception("Configuration") << "standard is temporarily disabled";
93  } else if (build == "cloneEngine") {
94  //buildFunction_ = mkfit::runBuildingTestPlexCloneEngine;
95  } else {
96  throw cms::Exception("Configuration")
97  << "Invalid value for parameter 'buildingRoutine' " << build << ", allowed are bestHit, standard, cloneEngine";
98  }
99 
100  // TODO: what to do when we have multiple instances of MkFitProducer in a job?
102 }
103 
105 
108 
109  desc.add("pixelHits", edm::InputTag("mkFitSiPixelHits"));
110  desc.add("stripHits", edm::InputTag("mkFitSiStripHits"));
111  desc.add("eventOfHits", edm::InputTag("mkFitEventOfHits"));
112  desc.add("seeds", edm::InputTag("mkFitSeedConverter"));
113  desc.add("clustersToSkip", edm::InputTag());
114  desc.add<std::string>("buildingRoutine", "cloneEngine")
115  ->setComment("Valid values are: 'bestHit', 'standard', 'cloneEngine'");
116  desc.add<edm::ESInputTag>("config")->setComment(
117  "ESProduct that has the mkFit configuration parameters for this iteration");
118  desc.add("seedCleaning", true)->setComment("Clean seeds within mkFit");
119  desc.add("removeDuplicates", true)->setComment("Run duplicate removal within mkFit");
120  desc.add("backwardFitInCMSSW", false)
121  ->setComment("Do backward fit (to innermost hit) in CMSSW (true) or mkFit (false)");
122  desc.addUntracked("mkFitSilent", true)->setComment("Allows to enables printouts from mkFit with 'False'");
123  desc.addUntracked("limitConcurrency", false)
124  ->setComment(
125  "Use tbb::task_arena to limit the internal concurrency to 1; useful only for timing studies when measuring "
126  "the module time");
127 
129  descCCC.add<double>("value");
130  desc.add("minGoodStripCharge", descCCC);
131 
132  descriptions.add("mkFitProducerDefault", desc);
133 }
134 
135 std::unique_ptr<mkfit::MkBuilderWrapper> MkFitProducer::beginStream(edm::StreamID iID) const {
136  return std::make_unique<mkfit::MkBuilderWrapper>(mkFitSilent_);
137 }
138 
140  const auto& pixelHits = iEvent.get(pixelHitsToken_);
141  const auto& stripHits = iEvent.get(stripHitsToken_);
142  const auto& eventOfHits = iEvent.get(eventOfHitsToken_);
143  const auto& seeds = iEvent.get(seedToken_);
144  if (seeds.seeds().empty()) {
146  return;
147  }
148  // This producer does not strictly speaking need the MkFitGeometry,
149  // but the ESProducer sets global variables (yes, that "feature"
150  // should be removed), so getting the MkFitGeometry makes it
151  // sure that the ESProducer is called even if the input/output
152  // converters
153  const auto& mkFitGeom = iSetup.getData(mkFitGeomToken_);
154  const auto& mkFitIterConfig = iSetup.getData(mkFitIterConfigToken_);
155 
156  const std::vector<bool>* pixelMaskPtr = nullptr;
157  std::vector<bool> pixelMask;
158  std::vector<bool> stripMask(stripHits.hits().size(), false);
159  if (not pixelMaskToken_.isUninitialized()) {
160  if (not pixelHits.hits().empty()) {
161  const auto& pixelContainerMask = iEvent.get(pixelMaskToken_);
162  pixelMask.resize(pixelContainerMask.size(), false);
163  if UNLIKELY (pixelContainerMask.refProd().id() != pixelHits.clustersID()) {
164  throw cms::Exception("LogicError") << "MkFitHitWrapper has pixel cluster ID " << pixelHits.clustersID()
165  << " but pixel cluster mask has " << pixelContainerMask.refProd().id();
166  }
167  pixelContainerMask.copyMaskTo(pixelMask);
168  pixelMaskPtr = &pixelMask;
169  }
170 
171  if (not stripHits.hits().empty()) {
172  const auto& stripContainerMask = iEvent.get(stripMaskToken_);
173  if UNLIKELY (stripContainerMask.refProd().id() != stripHits.clustersID()) {
174  throw cms::Exception("LogicError") << "MkFitHitWrapper has strip cluster ID " << stripHits.clustersID()
175  << " but strip cluster mask has " << stripContainerMask.refProd().id();
176  }
177  stripContainerMask.copyMaskTo(stripMask);
178  }
179  } else {
180  if (mkFitGeom.isPhase1())
182  }
183 
184  // seeds need to be mutable because of the possible cleaning
185  auto seeds_mutable = seeds.seeds();
187 
188  auto lambda = [&]() {
189  mkfit::run_OneIteration(mkFitGeom.trackerInfo(),
190  mkFitIterConfig,
191  eventOfHits.get(),
192  {pixelMaskPtr, &stripMask},
193  streamCache(iID)->get(),
194  seeds_mutable,
195  tracks,
199  };
200 
201  if (limitConcurrency_) {
202  tbb::task_arena arena(1);
203  arena.execute(std::move(lambda));
204  } else {
205  tbb::this_task_arena::isolate(std::move(lambda));
206  }
207 
209 }
210 
211 void MkFitProducer::stripClusterChargeCut(const std::vector<float>& stripClusterCharge, std::vector<bool>& mask) const {
212  if (mask.size() != stripClusterCharge.size()) {
213  throw cms::Exception("LogicError") << "Mask size (" << mask.size() << ") inconsistent with number of hits ("
214  << stripClusterCharge.size() << ")";
215  }
216  for (int i = 0, end = stripClusterCharge.size(); i < end; ++i) {
217  // mask == true means skip the cluster
218  mask[i] = mask[i] || (stripClusterCharge[i] <= minGoodStripCharge_);
219  }
220 }
221 
const bool seedCleaning_
edm::EDGetTokenT< edm::ContainerMask< edmNew::DetSetVector< SiStripCluster > > > stripMaskToken_
~MkFitProducer() override
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
edm::EDGetTokenT< edm::ContainerMask< edmNew::DetSetVector< SiPixelCluster > > > pixelMaskToken_
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
const bool backwardFitInCMSSW_
const bool removeDuplicates_
const edm::EDPutTokenT< MkFitOutputWrapper > putToken_
constexpr bool isUninitialized() const noexcept
Definition: EDGetToken.h:98
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const edm::EDGetTokenT< std::vector< float > > stripClusterChargeToken_
const bool limitConcurrency_
void stripClusterChargeCut(const std::vector< float > &stripClusterCharge, std::vector< bool > &mask) const
void run_OneIteration(const TrackerInfo &trackerInfo, const IterationConfig &itconf, const EventOfHits &eoh, const std::vector< const std::vector< bool > *> &hit_masks, MkBuilder &builder, TrackVec &seeds, TrackVec &out_tracks, bool do_seed_clean, bool do_backward_fit, bool do_remove_duplicates)
Definition: runFunctions.cc:27
const edm::EDGetTokenT< MkFitHitWrapper > stripHitsToken_
int iEvent
Definition: GenABIO.cc:224
std::unique_ptr< mkfit::MkBuilderWrapper > beginStream(edm::StreamID) const override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ParameterDescriptionBase * add(U const &iLabel, T const &value)
const float minGoodStripCharge_
std::vector< Track > TrackVec
void add(std::string const &label, ParameterSetDescription const &psetDescription)
const bool mkFitSilent_
const edm::EDGetTokenT< MkFitHitWrapper > pixelHitsToken_
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
const edm::ESGetToken< mkfit::IterationConfig, TrackerRecoGeometryRecord > mkFitIterConfigToken_
const edm::ESGetToken< MkFitGeometry, TrackerRecoGeometryRecord > mkFitGeomToken_
const edm::EDGetTokenT< MkFitEventOfHits > eventOfHitsToken_
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
#define UNLIKELY(x)
Definition: Likely.h:21
MkFitProducer(edm::ParameterSet const &iConfig)
const edm::EDGetTokenT< MkFitSeedWrapper > seedToken_
def move(src, dest)
Definition: eostools.py:511