CMS 3D CMS Logo

HIProtoTrackFilterProducer.cc
Go to the documentation of this file.
3 
6 
14 
20 
24 
26 public:
27  explicit HIProtoTrackFilterProducer(const edm::ParameterSet& iConfig);
28  ~HIProtoTrackFilterProducer() override;
29 
30  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
31 
32 private:
33  void produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const override;
34 
37  double theTIPMax;
40 };
41 
43  theBeamSpotToken(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"))),
44  theSiPixelRecHitsToken(consumes<SiPixelRecHitCollection>(iConfig.getParameter<edm::InputTag>("siPixelRecHits"))),
45  theTIPMax(iConfig.getParameter<double>("tipMax")),
46  theChi2Max(iConfig.getParameter<double>("chi2")),
47  thePtMin(iConfig.getParameter<double>("ptMin")),
48  doVariablePtMin(iConfig.getParameter<bool>("doVariablePtMin"))
49 {
50  produces<PixelTrackFilter>();
51 }
52 
55 
56  desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
57  desc.add<edm::InputTag>("siPixelRecHits", edm::InputTag("siPixelRecHits"));
58  desc.add<double>("ptMin", 1.0);
59  desc.add<double>("tipMax", 1.0);
60  desc.add<double>("chi2", 1000);
61  desc.add<bool>("doVariablePtMin", true);
62 
63  descriptions.add("hiProtoTrackFilter", desc);
64 }
65 
67 
69  // Get the beam spot
71  iEvent.getByToken( theBeamSpotToken, bsHandle);
72  const reco::BeamSpot *beamSpot = bsHandle.product();
73 
74  if(beamSpot) {
75  edm::LogInfo("HeavyIonVertexing")
76  << "[HIProtoTrackFilterProducer] Proto track selection based on beamspot"
77  << "\n (x,y,z) = (" << beamSpot->x0() << "," << beamSpot->y0() << "," << beamSpot->z0() << ")";
78  } else {
81  edm::LogError("HeavyIonVertexing") // this can be made a warning when operator() is fixed
82  << "No beamspot found with tag '" << labels.module << "'";
83  }
84 
85  // Estimate multiplicity
87  iEvent.getByToken(theSiPixelRecHitsToken, recHitColl);
88 
90  iSetup.get<TrackerTopologyRcd>().get(httopo);
91 
92  std::vector<const TrackingRecHit*> theChosenHits;
93  edmNew::copyDetSetRange(*recHitColl,theChosenHits, httopo->pxbDetIdLayerComparator(1));
94  float estMult = theChosenHits.size();
95 
96  double variablePtMin=thePtMin;
97  if(doVariablePtMin) {
98  // parameterize ptMin such that a roughly constant number of selected prototracks passed are to vertexing
99  float varPtCutoff = 1500; //cutoff for variable ptMin
100  if(estMult < varPtCutoff) {
101  variablePtMin = 0.075;
102  if(estMult > 0) variablePtMin = (13. - (varPtCutoff/estMult) )/12.;
103  if(variablePtMin<0.075) variablePtMin = 0.075; // don't lower the cut past 75 MeV
104  }
105  LogTrace("heavyIonHLTVertexing")<<" [HIProtoTrackFilterProducer: variablePtMin: " << variablePtMin << "]";
106  }
107 
108 
109  auto impl = std::make_unique<HIProtoTrackFilter>(beamSpot, theTIPMax, theChi2Max, variablePtMin);
110  auto prod = std::make_unique<PixelTrackFilter>(std::move(impl));
111  iEvent.put(std::move(prod));
112 }
113 
double z0() const
z coordinate
Definition: BeamSpot.h:68
void copyDetSetRange(DSTV const &dstv, std::vector< T const * > &v, std::pair< A, B > const &sel)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
char const * module
Definition: ProductLabels.h:5
edm::EDGetTokenT< SiPixelRecHitCollection > theSiPixelRecHitsToken
#define LogTrace(id)
T const * product() const
Definition: Handle.h:74
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
fixed size matrix
HLT enums.
T get() const
Definition: EventSetup.h:71
double y0() const
y coordinate
Definition: BeamSpot.h:66
edm::EDGetTokenT< reco::BeamSpot > theBeamSpotToken
std::pair< DetId, SameLayerComparator > pxbDetIdLayerComparator(uint32_t layer) const
HIProtoTrackFilterProducer(const edm::ParameterSet &iConfig)
void produce(edm::StreamID, edm::Event &iEvent, const edm::EventSetup &iSetup) const override
def move(src, dest)
Definition: eostools.py:511
double x0() const
x coordinate
Definition: BeamSpot.h:64