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 sumw = 0.;
30 
31  for (const auto& trk : vtx.originalTracks()) {
32  const double time = trk.timeExt();
33  const double err = trk.dtErrorExt();
35  continue; // tracks with no time information, as implemented in TransientTrackBuilder.cc l.17
36  const double inverr = err > 0. ? 1.0 / err : 0.;
37  const double w = inverr * inverr;
38  sumwt += w * time;
39  sumw += w;
40  }
41 
42  if (sumw > 0) {
43  vtxTime = sumwt / sumw;
44  vtxTimeError = 1 / sqrt(sumw);
45  return true;
46  }
47 
48  vtxTime = 0;
49  vtxTimeError = 1.;
50  return false;
51 }
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:23
static constexpr float defaultInvalidTrackTimeReso
void setEvent(edm::Event &iEvent, edm::EventSetup const &iSetup) override