CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
GlobalTrackingRegionProducerFromBeamSpot.h
Go to the documentation of this file.
1 #ifndef RecoTracker_TkTrackingRegions_GlobalTrackingRegionProducerFromBeamSpot_H
2 #define RecoTracker_TkTrackingRegions_GlobalTrackingRegionProducerFromBeamSpot_H
3 
17 
19 public:
23  edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
24  thePtMin = regionPSet.getParameter<double>("ptMin");
25  theOriginRadius = regionPSet.getParameter<double>("originRadius");
26  if (!regionPSet.existsAs<double>("nSigmaZ") && !regionPSet.existsAs<double>("originHalfLength")) {
27  throw cms::Exception("Configuration") << "GlobalTrackingRegionProducerFromBeamSpot: at least one of nSigmaZ, "
28  "originHalfLength must be present in the cfg.\n";
29  }
30  theNSigmaZ = (regionPSet.existsAs<double>("nSigmaZ") ? regionPSet.getParameter<double>("nSigmaZ") : 0.0);
32  (regionPSet.existsAs<double>("originHalfLength") ? regionPSet.getParameter<double>("originHalfLength") : 0.0);
33  token_beamSpot = iC.consumes<reco::BeamSpot>(regionPSet.getParameter<edm::InputTag>("beamSpot"));
34  thePrecise = regionPSet.getParameter<bool>("precise");
35  theUseMS =
36  (regionPSet.existsAs<bool>("useMultipleScattering") ? regionPSet.getParameter<bool>("useMultipleScattering")
37  : false);
38  if (theUseMS) {
40  }
41  }
42 
43  ~GlobalTrackingRegionProducerFromBeamSpot() override = default;
44 
45  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
46  {
48 
49  desc.add<bool>("precise", true);
50  desc.add<bool>("useMultipleScattering", false);
51  desc.add<double>("nSigmaZ", 4.0);
52  desc.add<double>("originHalfLength", 0.0); // this is the default in constructor
53  desc.add<double>("originRadius", 0.2);
54  desc.add<double>("ptMin", 0.9);
55  desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
56 
57  // Only for backwards-compatibility
59  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
60 
61  descriptions.add("globalTrackingRegionFromBeamSpot", descRegion);
62  }
63 
64  {
66 
67  desc.add<bool>("precise", true);
68  desc.add<bool>("useMultipleScattering", false);
69  desc.add<double>("nSigmaZ", 0.0); // this is the default in constructor
70  desc.add<double>("originHalfLength", 21.2);
71  desc.add<double>("originRadius", 0.2);
72  desc.add<double>("ptMin", 0.9);
73  desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
74 
75  // Only for backwards-compatibility
77  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
78 
79  descriptions.add("globalTrackingRegionFromBeamSpotFixedZ", descRegion);
80  }
81  }
82 
83  std::vector<std::unique_ptr<TrackingRegion> > regions(const edm::Event& ev,
84  const edm::EventSetup& es) const override {
85  std::vector<std::unique_ptr<TrackingRegion> > result;
87  ev.getByToken(token_beamSpot, bsHandle);
88  if (bsHandle.isValid()) {
89  const reco::BeamSpot& bs = *bsHandle;
90 
91  GlobalPoint origin(bs.x0(), bs.y0(), bs.z0());
92 
93  const MultipleScatteringParametrisationMaker* msmaker = nullptr;
94  if (theUseMS) {
95  msmaker = &es.getData(token_msmaker);
96  }
97 
98  result.push_back(std::make_unique<GlobalTrackingRegion>(thePtMin,
99  origin,
102  thePrecise,
103  theUseMS,
104  msmaker));
105  }
106  return result;
107  }
108 
109 private:
110  double thePtMin;
113  double theNSigmaZ;
117  bool theUseMS;
118 };
119 
120 #endif
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
double z0() const
z coordinate
Definition: BeamSpot.h:65
tuple cfg
Definition: looper.py:296
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:171
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
~GlobalTrackingRegionProducerFromBeamSpot() override=default
GlobalTrackingRegionProducerFromBeamSpot(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
tuple result
Definition: mps_fire.py:311
bool getData(T &iHolder) const
Definition: EventSetup.h:122
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool isValid() const
Definition: HandleBase.h:70
GlobalTrackingRegionProducerFromBeamSpot(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
double sigmaZ() const
sigma z
Definition: BeamSpot.h:76
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > token_msmaker
double y0() const
y coordinate
Definition: BeamSpot.h:63
std::vector< std::unique_ptr< TrackingRegion > > regions(const edm::Event &ev, const edm::EventSetup &es) const override
double x0() const
x coordinate
Definition: BeamSpot.h:61