CMS 3D CMS Logo

VertexTimeAlgorithmLegacy4D.cc
Go to the documentation of this file.
6 
7 #ifdef PVTX_DEBUG
8 #define LOG edm::LogPrint("VertexTimeAlgorithmLegacy4D")
9 #else
10 #define LOG LogDebug("VertexTimeAlgorithmLegacy4D")
11 #endif
12 
14  : VertexTimeAlgorithmBase(iConfig, iCC) {}
15 
18 }
19 
21 
22 bool VertexTimeAlgorithmLegacy4D::vertexTime(float& vtxTime, float& vtxTimeError, const TransientVertex& vtx) const {
23  const auto num_track = vtx.originalTracks().size();
24  if (num_track == 0) {
25  return false;
26  }
27 
28  double sumwt = 0.;
29  double sumwt2 = 0.;
30  double sumw = 0.;
31  double vartime = 0.;
32 
33  for (const auto& trk : vtx.originalTracks()) {
34  const double time = trk.timeExt();
35  const double err = trk.dtErrorExt();
37  continue; // tracks with no time information, as implemented in TransientTrackBuilder.cc l.17
38  const double inverr = err > 0. ? 1.0 / err : 0.;
39  const double w = inverr * inverr;
40  sumwt += w * time;
41  sumwt2 += w * time * time;
42  sumw += w;
43  }
44 
45  if (sumw > 0) {
46  double sumsq = sumwt2 - sumwt * sumwt / sumw;
47  double chisq = num_track > 1 ? sumsq / double(num_track - 1) : sumsq / double(num_track);
48  vartime = chisq / sumw;
49 
50  vtxTime = sumwt / sumw;
51  vtxTimeError = sqrt(vartime);
52  return true;
53  }
54 
55  vtxTime = 0;
56  vtxTimeError = 1.;
57  return false;
58 }
static void fillPSetDescription(edm::ParameterSetDescription &iDesc)
T w() const
VertexTimeAlgorithmLegacy4D(const edm::ParameterSet &conf, edm::ConsumesCollector &iC)
bool vertexTime(float &vtxTime, float &vtxTimeError, TransientVertex const &vtx) const override
static void fillPSetDescription(edm::ParameterSetDescription &iDesc)
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
static constexpr float defaultInvalidTrackTimeReso
void setEvent(edm::Event &iEvent, edm::EventSetup const &iSetup) override