CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Attributes
GlobalTrackingRegionWithVerticesProducer Class Reference

#include <GlobalTrackingRegionWithVerticesProducer.h>

Inheritance diagram for GlobalTrackingRegionWithVerticesProducer:
TrackingRegionProducer

Public Member Functions

 GlobalTrackingRegionWithVerticesProducer (const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
 
std::vector< std::unique_ptr< TrackingRegion > > regions (const edm::Event &ev, const edm::EventSetup &) const override
 
 ~GlobalTrackingRegionWithVerticesProducer () override
 
- Public Member Functions inherited from TrackingRegionProducer
virtual ~TrackingRegionProducer ()
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 

Private Attributes

edm::InputTag theBeamSpotTag
 
double theFixedError
 
bool theHalfLengthScaling
 
int theMaxNVertices
 
double theMaxPtMin
 
double theMinHalfLength
 
double theMinOriginR
 
double theNSigmaZ
 
double theOriginRadius
 
bool theOriginRScaling
 
bool thePrecise
 
double thePtMin
 
bool thePtMinScaling
 
double theScalingEnd
 
double theScalingStart
 
double theSigmaZVertex
 
bool theUseFakeVertices
 
bool theUseFixedError
 
bool theUseFoundVertices
 
bool theUseMS
 
edm::EDGetTokenT< reco::BeamSpottoken_beamSpot
 
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > token_pc
 
edm::EDGetTokenT< reco::VertexCollectiontoken_vertex
 

Detailed Description

Definition at line 18 of file GlobalTrackingRegionWithVerticesProducer.h.

Constructor & Destructor Documentation

GlobalTrackingRegionWithVerticesProducer::GlobalTrackingRegionWithVerticesProducer ( const edm::ParameterSet cfg,
edm::ConsumesCollector &&  iC 
)
inline

Definition at line 20 of file GlobalTrackingRegionWithVerticesProducer.h.

References edm::ParameterSet::getParameter(), HLT_FULL_cff::pixelClustersForScaling, theFixedError, theHalfLengthScaling, theMaxNVertices, theMaxPtMin, theMinHalfLength, theMinOriginR, theNSigmaZ, theOriginRadius, theOriginRScaling, thePrecise, thePtMin, thePtMinScaling, theScalingEnd, theScalingStart, theSigmaZVertex, theUseFakeVertices, theUseFixedError, theUseFoundVertices, theUseMS, token_beamSpot, token_pc, and token_vertex.

20  {
21  edm::ParameterSet regionPSet = cfg.getParameter<edm::ParameterSet>("RegionPSet");
22 
23  thePtMin = regionPSet.getParameter<double>("ptMin");
24  theOriginRadius = regionPSet.getParameter<double>("originRadius");
25  theNSigmaZ = regionPSet.getParameter<double>("nSigmaZ");
26  token_beamSpot = iC.consumes<reco::BeamSpot>(regionPSet.getParameter<edm::InputTag>("beamSpot"));
27  thePrecise = regionPSet.getParameter<bool>("precise");
28  theUseMS = regionPSet.getParameter<bool>("useMultipleScattering");
29 
30  theSigmaZVertex = regionPSet.getParameter<double>("sigmaZVertex");
31  theFixedError = regionPSet.getParameter<double>("fixedError");
32 
33  theMaxNVertices = regionPSet.getParameter<int>("maxNVertices");
34 
35  theUseFoundVertices = regionPSet.getParameter<bool>("useFoundVertices");
36  theUseFakeVertices = regionPSet.getParameter<bool>("useFakeVertices");
37  theUseFixedError = regionPSet.getParameter<bool>("useFixedError");
38  token_vertex = iC.consumes<reco::VertexCollection>(regionPSet.getParameter<edm::InputTag>("VertexCollection"));
39 
40  //information for Heavy ion region scaling
41  theOriginRScaling = regionPSet.getParameter<bool>("originRScaling4BigEvts");
42  thePtMinScaling = regionPSet.getParameter<bool>("ptMinScaling4BigEvts");
43  theHalfLengthScaling = regionPSet.getParameter<bool>("halfLengthScaling4BigEvts");
44  theMinOriginR = regionPSet.getParameter<double>("minOriginR");
45  theMaxPtMin = regionPSet.getParameter<double>("maxPtMin");
46  theMinHalfLength = regionPSet.getParameter<double>("minHalfLength");
47  theScalingStart = regionPSet.getParameter<double>("scalingStartNPix");
48  theScalingEnd = regionPSet.getParameter<double>("scalingEndNPix");
49  edm::InputTag pixelClustersForScaling = regionPSet.getParameter<edm::InputTag>("pixelClustersForScaling");
52  }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::VertexCollection > token_vertex
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > token_pc
GlobalTrackingRegionWithVerticesProducer::~GlobalTrackingRegionWithVerticesProducer ( )
inlineoverride

Definition at line 54 of file GlobalTrackingRegionWithVerticesProducer.h.

54 {}

Member Function Documentation

static void GlobalTrackingRegionWithVerticesProducer::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
inlinestatic

Definition at line 56 of file GlobalTrackingRegionWithVerticesProducer.h.

References edm::ConfigurationDescriptions::add(), edm::ParameterSetDescription::add(), and HLT_2018_cff::InputTag.

56  {
58 
59  desc.add<bool>("precise", true);
60  desc.add<bool>("useMultipleScattering", false);
61  desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
62  desc.add<bool>("useFixedError", true);
63  desc.add<double>("originRadius", 0.2);
64  desc.add<double>("sigmaZVertex", 3.0);
65  desc.add<double>("fixedError", 0.2);
66  desc.add<edm::InputTag>("VertexCollection", edm::InputTag("firstStepPrimaryVertices"));
67  desc.add<double>("ptMin", 0.9);
68  desc.add<bool>("useFoundVertices", true);
69  desc.add<bool>("useFakeVertices", false);
70  desc.add<int>("maxNVertices", -1)->setComment("-1 for all vertices");
71  desc.add<double>("nSigmaZ", 4.0);
72  desc.add<edm::InputTag>("pixelClustersForScaling", edm::InputTag("siPixelClusters"));
73  desc.add<bool>("originRScaling4BigEvts", false);
74  desc.add<bool>("ptMinScaling4BigEvts", false);
75  desc.add<bool>("halfLengthScaling4BigEvts", false);
76  desc.add<double>("minOriginR", 0);
77  desc.add<double>("maxPtMin", 1000);
78  desc.add<double>("minHalfLength", 0);
79  desc.add<double>("scalingStartNPix", 0.0);
80  desc.add<double>("scalingEndNPix", 1.0);
81 
82  // Only for backwards-compatibility
84  descRegion.add<edm::ParameterSetDescription>("RegionPSet", desc);
85 
86  descriptions.add("globalTrackingRegionWithVertices", descRegion);
87  }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
std::vector<std::unique_ptr<TrackingRegion> > GlobalTrackingRegionWithVerticesProducer::regions ( const edm::Event ev,
const edm::EventSetup  
) const
inlineoverridevirtual

Implements TrackingRegionProducer.

Definition at line 89 of file GlobalTrackingRegionWithVerticesProducer.h.

References edmNew::DetSetVector< T >::dataSize(), Exception, edm::Event::getByToken(), input, edm::HandleBase::isValid(), mps_fire::result, reco::BeamSpot::sigmaZ(), theFixedError, theHalfLengthScaling, theMaxNVertices, theMaxPtMin, theMinHalfLength, theMinOriginR, theNSigmaZ, theOriginRadius, theOriginRScaling, thePrecise, thePtMin, thePtMinScaling, theScalingEnd, theScalingStart, theSigmaZVertex, theUseFakeVertices, theUseFixedError, theUseFoundVertices, theUseMS, token_beamSpot, token_pc, token_vertex, spclusmultinvestigator_cfi::vertexCollection, reco::BeamSpot::x0(), reco::BeamSpot::y0(), and reco::BeamSpot::z0().

89  {
90  std::vector<std::unique_ptr<TrackingRegion> > result;
91 
92  GlobalPoint theOrigin;
94  ev.getByToken(token_beamSpot, bsHandle);
95  double bsSigmaZ;
96  if (bsHandle.isValid()) {
97  const reco::BeamSpot& bs = *bsHandle;
98  bsSigmaZ = theNSigmaZ * bs.sigmaZ();
99  theOrigin = GlobalPoint(bs.x0(), bs.y0(), bs.z0());
100  } else {
101  throw cms::Exception("Seeding") << "ERROR: input beamSpot is not valid in GlobalTrackingRegionWithVertices";
102  }
103 
104  if (theUseFoundVertices) {
106  ev.getByToken(token_vertex, vertexCollection);
107 
110  if (doScaling)
111  ev.getByToken(token_pc, pixelClusterDSV);
112 
113  for (reco::VertexCollection::const_iterator iV = vertexCollection->begin(); iV != vertexCollection->end(); iV++) {
114  if (!iV->isValid())
115  continue;
116  if (iV->isFake() && !(theUseFakeVertices && theUseFixedError))
117  continue;
118  GlobalPoint theOrigin_ = GlobalPoint(iV->x(), iV->y(), iV->z());
119 
120  //scaling origin radius, half length, min pt for high-occupancy HI events to keep timing reasonable
121  if (doScaling) {
122  //Use the unscaled radius unless one of the two conditions below is met
123  double scaledOriginRadius = theOriginRadius;
124  double scaledHalfLength = theFixedError;
125  double scaledPtMin = thePtMin;
126 
127  //calculate nPixels (adapted from TkSeedGenerator/src/ClusterChecker.cc)
128  double nPix = 0;
129 
130  const edmNew::DetSetVector<SiPixelCluster>& input = *pixelClusterDSV;
131  nPix = input.dataSize();
132 
133  //first condition is for high occupancy, second makes sure we won't divide by zero or a negative number
134  if ((nPix > theScalingEnd) || ((theScalingEnd - theScalingStart) <= 0)) {
135  if (theOriginRScaling)
136  scaledOriginRadius = theMinOriginR; // sets parameters to minimum value from PSet
138  scaledHalfLength = theMinHalfLength;
139  if (thePtMinScaling)
140  scaledPtMin = theMaxPtMin;
141  }
142  //second condition - scale radius linearly by Npix in the region from ScalingStart to ScalingEnd
143  else if ((nPix <= theScalingEnd) && (nPix > theScalingStart)) {
144  float slopeFactor = (nPix - theScalingStart) / (theScalingEnd - theScalingStart);
145  if (theOriginRScaling)
146  scaledOriginRadius = theOriginRadius - (theOriginRadius - theMinOriginR) * slopeFactor;
148  scaledHalfLength = theFixedError - (theFixedError - theMinHalfLength) * slopeFactor;
149  if (thePtMinScaling)
150  scaledPtMin = thePtMin - (thePtMin - theMaxPtMin) * slopeFactor;
151  }
152  //if region has 0 size, return 'result' empty, otherwise make a tracking region
153  if (scaledOriginRadius != 0 && scaledHalfLength != 0) {
154  result.push_back(std::make_unique<GlobalTrackingRegion>(
155  scaledPtMin, theOrigin_, scaledOriginRadius, scaledHalfLength, thePrecise, theUseMS));
156  }
157  } //end of region scaling code, pp behavior below
158 
159  else {
160  double theOriginHalfLength_ = (theUseFixedError ? theFixedError : (iV->zError()) * theSigmaZVertex);
161  result.push_back(std::make_unique<GlobalTrackingRegion>(
162  thePtMin, theOrigin_, theOriginRadius, theOriginHalfLength_, thePrecise, theUseMS));
163  if (theMaxNVertices >= 0 && result.size() >= static_cast<unsigned>(theMaxNVertices))
164  break;
165  }
166  }
167 
168  if (result.empty() && !(theOriginRScaling || thePtMinScaling || theHalfLengthScaling)) {
169  result.push_back(std::make_unique<GlobalTrackingRegion>(
170  thePtMin, theOrigin, theOriginRadius, bsSigmaZ, thePrecise, theUseMS));
171  }
172  } else {
173  result.push_back(
174  std::make_unique<GlobalTrackingRegion>(thePtMin, theOrigin, theOriginRadius, bsSigmaZ, thePrecise, theUseMS));
175  }
176 
177  return result;
178  }
double z0() const
z coordinate
Definition: BeamSpot.h:65
size_type dataSize() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
edm::EDGetTokenT< reco::VertexCollection > token_vertex
edm::EDGetTokenT< edmNew::DetSetVector< SiPixelCluster > > token_pc
static std::string const input
Definition: EdmProvDump.cc:48
bool isValid() const
Definition: HandleBase.h:70
double sigmaZ() const
sigma z
Definition: BeamSpot.h:76
double y0() const
y coordinate
Definition: BeamSpot.h:63
double x0() const
x coordinate
Definition: BeamSpot.h:61

Member Data Documentation

edm::InputTag GlobalTrackingRegionWithVerticesProducer::theBeamSpotTag
private

Definition at line 184 of file GlobalTrackingRegionWithVerticesProducer.h.

double GlobalTrackingRegionWithVerticesProducer::theFixedError
private
bool GlobalTrackingRegionWithVerticesProducer::theHalfLengthScaling
private
int GlobalTrackingRegionWithVerticesProducer::theMaxNVertices
private
double GlobalTrackingRegionWithVerticesProducer::theMaxPtMin
private
double GlobalTrackingRegionWithVerticesProducer::theMinHalfLength
private
double GlobalTrackingRegionWithVerticesProducer::theMinOriginR
private
double GlobalTrackingRegionWithVerticesProducer::theNSigmaZ
private
double GlobalTrackingRegionWithVerticesProducer::theOriginRadius
private
bool GlobalTrackingRegionWithVerticesProducer::theOriginRScaling
private
bool GlobalTrackingRegionWithVerticesProducer::thePrecise
private
double GlobalTrackingRegionWithVerticesProducer::thePtMin
private
bool GlobalTrackingRegionWithVerticesProducer::thePtMinScaling
private
double GlobalTrackingRegionWithVerticesProducer::theScalingEnd
private
double GlobalTrackingRegionWithVerticesProducer::theScalingStart
private
double GlobalTrackingRegionWithVerticesProducer::theSigmaZVertex
private
bool GlobalTrackingRegionWithVerticesProducer::theUseFakeVertices
private
bool GlobalTrackingRegionWithVerticesProducer::theUseFixedError
private
bool GlobalTrackingRegionWithVerticesProducer::theUseFoundVertices
private
bool GlobalTrackingRegionWithVerticesProducer::theUseMS
private
edm::EDGetTokenT<reco::BeamSpot> GlobalTrackingRegionWithVerticesProducer::token_beamSpot
private
edm::EDGetTokenT<edmNew::DetSetVector<SiPixelCluster> > GlobalTrackingRegionWithVerticesProducer::token_pc
private
edm::EDGetTokenT<reco::VertexCollection> GlobalTrackingRegionWithVerticesProducer::token_vertex
private