CMS 3D CMS Logo

HIBestVertexProducer.cc
Go to the documentation of this file.
2 
8 
10 
11 #include <algorithm>
12 #include <iostream>
13 using namespace std;
14 using namespace edm;
15 
16 /*****************************************************************************/
18  : theConfig(ps),
19  theBeamSpotTag(consumes<reco::BeamSpot>(ps.getParameter<edm::InputTag>("beamSpotLabel"))),
20  theMedianVertexCollection(
21  consumes<reco::VertexCollection>(ps.getParameter<edm::InputTag>("medianVertexCollection"))),
22  theAdaptiveVertexCollection(
23  consumes<reco::VertexCollection>(ps.getParameter<edm::InputTag>("adaptiveVertexCollection"))),
24  theUseFinalAdaptiveVertexCollection(ps.getParameter<bool>("useFinalAdaptiveVertexCollection")) {
27  consumes<reco::VertexCollection>(ps.getParameter<edm::InputTag>("finalAdaptiveVertexCollection"));
28  }
29  produces<reco::VertexCollection>();
30 }
31 
32 /*****************************************************************************/
34 
35 /*****************************************************************************/
38  desc.add<InputTag>("beamSpotLabel", edm::InputTag("offlineBeamSpot"));
39  desc.add<InputTag>("adaptiveVertexCollection", edm::InputTag("hiBestAdaptiveVertex"));
40  desc.add<InputTag>("medianVertexCollection", edm::InputTag("hiPixelMedianVertex"));
41  desc.add<bool>("useFinalAdaptiveVertexCollection", false);
42  desc.add<InputTag>("finalAdaptiveVertexCollection", edm::InputTag("hiBestOfflinePrimaryVertex"));
43 
44  descriptions.add("hiSelectedPixelVertex", desc);
45 }
46 
47 /*****************************************************************************/
49 
50 /*****************************************************************************/
52  // 1. use best adaptive vertex preferentially
53  // 2. use median vertex in case where adaptive algorithm failed
54  // 3. use beamspot if netither vertexing method succeeds
55 
56  // New vertex collection
57  auto newVertexCollection = std::make_unique<reco::VertexCollection>();
58 
59  //** Get precise adaptive vertex **/
61  ev.getByToken(theAdaptiveVertexCollection, vc1);
62  const reco::VertexCollection* vertices1 = vc1.product();
63 
64  if (vertices1->empty())
65  LogError("HeavyIonVertexing") << "adaptive vertex collection is empty!" << endl;
66 
67  //** Final vertex collection if needed **//
68  const double maxZError = 3.0; //any vtx with error > this number is considered no good
69  bool hasFinalVertex = false;
72  ev.getByToken(theFinalAdaptiveVertexCollection, vc0);
73  const reco::VertexCollection* vertices0 = vc0.product();
74  if (vertices0->empty())
75  LogInfo("HeavyIonVertexing") << "final adaptive vertex collection is empty!" << endl;
76 
77  //if using final vertex and has a good vertex, use it
78  if (vertices0->begin()->zError() < maxZError) {
79  hasFinalVertex = true;
80  auto const& vertex0 = vertices0->front();
81  newVertexCollection->push_back(vertex0);
82  LogInfo("HeavyIonVertexing") << "adaptive vertex:\n vz = (" << vertex0.x() << ", " << vertex0.y() << ", "
83  << vertex0.z() << ")"
84  << "\n error = (" << vertex0.xError() << ", " << vertex0.yError() << ", "
85  << vertex0.zError() << ")" << endl;
86  }
87  }
88 
89  //otherwise use the pixel track adaptive vertex if it is good
90  if (!hasFinalVertex) {
91  if (vertices1->begin()->zError() < maxZError) {
92  reco::VertexCollection::const_iterator vertex1 = vertices1->begin();
93  newVertexCollection->push_back(*vertex1);
94 
95  LogInfo("HeavyIonVertexing") << "adaptive vertex:\n vz = (" << vertex1->x() << ", " << vertex1->y() << ", "
96  << vertex1->z() << ")"
97  << "\n error = (" << vertex1->xError() << ", " << vertex1->yError() << ", "
98  << vertex1->zError() << ")" << endl;
99  } else {
100  //** Get fast median vertex **/
102  ev.getByToken(theMedianVertexCollection, vc2);
103  const reco::VertexCollection* vertices2 = vc2.product();
104 
105  //** Get beam spot position and error **/
107  edm::Handle<reco::BeamSpot> beamSpotHandle;
108  ev.getByToken(theBeamSpotTag, beamSpotHandle);
109 
110  if (beamSpotHandle.isValid())
111  beamSpot = *beamSpotHandle;
112  else
113  LogError("HeavyIonVertexing") << "no beamspot found " << endl;
114 
115  if (!vertices2->empty()) {
116  reco::VertexCollection::const_iterator vertex2 = vertices2->begin();
118  err(0, 0) = pow(beamSpot.BeamWidthX(), 2);
119  err(1, 1) = pow(beamSpot.BeamWidthY(), 2);
120  err(2, 2) = pow(vertex2->zError(), 2);
121  reco::Vertex newVertex(reco::Vertex::Point(beamSpot.x0(), beamSpot.y0(), vertex2->z()), err, 0, 1, 1);
122  newVertexCollection->push_back(newVertex);
123 
124  LogInfo("HeavyIonVertexing") << "median vertex + beamspot: \n position = (" << newVertex.x() << ", "
125  << newVertex.y() << ", " << newVertex.z() << ")"
126  << "\n error = (" << newVertex.xError() << ", " << newVertex.yError() << ", "
127  << newVertex.zError() << ")" << endl;
128 
129  } else {
131  err(0, 0) = pow(beamSpot.BeamWidthX(), 2);
132  err(1, 1) = pow(beamSpot.BeamWidthY(), 2);
133  err(2, 2) = pow(beamSpot.sigmaZ(), 2);
134  reco::Vertex newVertex(beamSpot.position(), err, 0, 0, 1);
135  newVertexCollection->push_back(newVertex);
136 
137  LogInfo("HeavyIonVertexing") << "beam spot: \n position = (" << newVertex.x() << ", " << newVertex.y() << ", "
138  << newVertex.z() << ")"
139  << "\n error = (" << newVertex.xError() << ", " << newVertex.yError() << ", "
140  << newVertex.zError() << ")" << endl;
141  }
142  }
143  }
144 
145  // put new vertex collection into event
146  ev.put(std::move(newVertexCollection));
147 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
T const * product() const
Definition: Handle.h:70
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
constexpr int pow(int x)
Definition: conifer.h:24
edm::EDGetTokenT< reco::VertexCollection > theAdaptiveVertexCollection
Log< level::Error, false > LogError
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
edm::EDGetTokenT< reco::BeamSpot > theBeamSpotTag
void produce(edm::Event &ev, const edm::EventSetup &es) override
edm::EDGetTokenT< reco::VertexCollection > theFinalAdaptiveVertexCollection
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Log< level::Info, false > LogInfo
HIBestVertexProducer(const edm::ParameterSet &ps)
edm::EDGetTokenT< reco::VertexCollection > theMedianVertexCollection
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool isValid() const
Definition: HandleBase.h:70
fixed size matrix
HLT enums.
def move(src, dest)
Definition: eostools.py:511