CMS 3D CMS Logo

BPHDecayVertex.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author Paolo Ronchese INFN Padova
5  *
6  */
7 
8 //-----------------------
9 // This Class' Header --
10 //-----------------------
12 
13 //-------------------------------
14 // Collaborating Class Headers --
15 //-------------------------------
30 
31 //---------------
32 // C++ Headers --
33 //---------------
34 using namespace std;
35 
36 //-------------------
37 // Initializations --
38 //-------------------
39 
40 //----------------
41 // Constructors --
42 //----------------
43 BPHDecayVertex::BPHDecayVertex(const BPHEventSetupWrapper* es, int daugNum, int compNum)
44  : BPHDecayMomentum(daugNum, compNum),
45  evSetup(new BPHEventSetupWrapper(es)),
46  oldTracks(true),
47  oldTTracks(true),
48  oldVertex(true),
49  validTks(false),
50  savedFitter(nullptr),
51  savedBS(nullptr),
52  savedPP(nullptr),
53  savedPE(nullptr) {}
54 
56  : evSetup(new BPHEventSetupWrapper(es)),
57  oldTracks(true),
58  oldTTracks(true),
59  oldVertex(true),
60  validTks(false),
61  savedFitter(nullptr),
62  savedBS(nullptr),
63  savedPP(nullptr),
64  savedPE(nullptr) {
65  map<const reco::Candidate*, const reco::Candidate*> iMap;
66  const vector<const reco::Candidate*>& daug = daughters();
67  const vector<Component>& list = ptr->componentList();
68  int i;
69  int n = daug.size();
70  for (i = 0; i < n; ++i) {
71  const reco::Candidate* cand = daug[i];
72  iMap[originalReco(cand)] = cand;
73  }
74  for (i = 0; i < n; ++i) {
75  const Component& c = list[i];
76  searchMap[iMap[c.cand]] = c.searchList;
77  }
78  const vector<BPHRecoConstCandPtr>& dComp = daughComp();
79  int j;
80  int m = dComp.size();
81  for (j = 0; j < m; ++j) {
82  const map<const reco::Candidate*, string>& dMap = dComp[j]->searchMap;
83  searchMap.insert(dMap.begin(), dMap.end());
84  }
85 }
86 
87 //--------------
88 // Destructor --
89 //--------------
91 
92 //--------------
93 // Operations --
94 //--------------
96  if (oldTracks)
97  fTracks();
98  return validTks;
99 }
100 
102  vertex();
103  return validTks && fittedVertex.isValid();
104 }
105 
107  const reco::BeamSpot* bs,
108  const GlobalPoint* priorPos,
109  const GlobalError* priorError) const {
110  if ((fitter == nullptr) && (bs == nullptr) && (priorPos == nullptr) && (priorError == nullptr)) {
111  fitter = savedFitter;
112  bs = savedBS;
113  priorPos = savedPP;
114  priorError = savedPE;
115  }
116  if (oldVertex || (fitter != savedFitter) || (bs != savedBS) || (priorPos != savedPP) || (priorError != savedPE)) {
117  if (fitter != nullptr) {
118  fitVertex(fitter, bs, priorPos, priorError);
119  } else {
120  KalmanVertexFitter kvf(true);
121  fitVertex(&kvf, bs, priorPos, priorError);
122  }
123  }
124  return fittedVertex;
125 }
126 
127 const vector<const reco::Track*>& BPHDecayVertex::tracks() const {
128  if (oldTracks)
129  fTracks();
130  return rTracks;
131 }
132 
134  if (oldTracks)
135  fTracks();
136  map<const reco::Candidate*, const reco::Track*>::const_iterator iter = tkMap.find(cand);
137  map<const reco::Candidate*, const reco::Track*>::const_iterator iend = tkMap.end();
138  return (iter != iend ? iter->second : nullptr);
139 }
140 
142  if (oldTracks)
143  fTracks();
144  map<const reco::Candidate*, char>::const_iterator iter = tmMap.find(cand);
145  map<const reco::Candidate*, char>::const_iterator iend = tmMap.end();
146  return (iter != iend ? iter->second : '.');
147 }
148 
149 const vector<reco::TransientTrack>& BPHDecayVertex::transientTracks() const {
150  if (oldTTracks)
151  fTTracks();
152  return trTracks;
153 }
154 
156  if (oldTTracks)
157  fTTracks();
158  map<const reco::Candidate*, reco::TransientTrack*>::const_iterator iter = ttMap.find(cand);
159  map<const reco::Candidate*, reco::TransientTrack*>::const_iterator iend = ttMap.end();
160  return (iter != iend ? iter->second : nullptr);
161 }
162 
165 
167  static const string dum = "";
168  map<const reco::Candidate*, string>::const_iterator iter = searchMap.find(cand);
169  if (iter != searchMap.end())
170  return iter->second;
171  return dum;
172 }
173 
174 void BPHDecayVertex::addV(const string& name, const reco::Candidate* daug, const string& searchList, double mass) {
175  addP(name, daug, mass);
176  searchMap[daughters().back()] = searchList;
177  return;
178 }
179 
180 void BPHDecayVertex::addV(const string& name, const BPHRecoConstCandPtr& comp) {
181  addP(name, comp);
182  const map<const reco::Candidate*, string>& dMap = comp->searchMap;
183  searchMap.insert(dMap.begin(), dMap.end());
184  return;
185 }
186 
189  oldTracks = oldVertex = true;
190  validTks = false;
191  return;
192 }
193 
195  oldTTracks = true;
196  rTracks.clear();
197  tkMap.clear();
198  tmMap.clear();
199  const vector<const reco::Candidate*>& dL = daughFull();
200  int n = dL.size();
201  trTracks.reserve(n);
202  validTks = true;
203  while (n--) {
204  const reco::Candidate* rp = dL[n];
205  tkMap[rp] = nullptr;
206  tmMap[rp] = '.';
207  if (!rp->charge())
208  continue;
209  const char* searchList = "cfhp";
210  char usedMode;
211  map<const reco::Candidate*, string>::const_iterator iter = searchMap.find(rp);
212  if (iter != searchMap.end())
213  searchList = iter->second.c_str();
214  const reco::Track* tp = tkMap[rp] = BPHTrackReference::getTrack(*originalReco(rp), searchList, &usedMode);
215  if (tp == nullptr) {
216  edm::LogPrint("DataNotFound") << "BPHDecayVertex::tTracks: "
217  << "no track for reco::(PF)Candidate";
218  validTks = false;
219  continue;
220  }
221  rTracks.push_back(tp);
222  tmMap[rp] = usedMode;
223  }
224  oldTracks = false;
225  return;
226 }
227 
229  if (oldTracks)
230  fTracks();
231  trTracks.clear();
234  const edm::EventSetup* ep = evSetup->get();
236  token->get(*ep, ttB);
237  ttMap.clear();
238  const vector<const reco::Candidate*>& dL = daughFull();
239  int n = dL.size();
240  trTracks.reserve(n);
241  while (n--) {
242  const reco::Candidate* rp = dL[n];
243  ttMap[rp] = nullptr;
244  map<const reco::Candidate*, const reco::Track*>::const_iterator iter = tkMap.find(rp);
245  const reco::Track* tp = iter->second;
246  if (tp == nullptr)
247  continue;
248  trTracks.push_back(ttB->build(tp));
249  ttMap[rp] = &trTracks.back();
250  }
251  oldTTracks = false;
252 }
253 
255  const reco::BeamSpot* bs,
256  const GlobalPoint* priorPos,
257  const GlobalError* priorError) const {
258  oldVertex = false;
259  savedFitter = fitter;
260  savedBS = bs;
261  savedPP = priorPos;
262  savedPE = priorError;
263  if (oldTTracks)
264  fTTracks();
265  if (trTracks.size() < 2)
266  return;
267  try {
268  if (bs == nullptr) {
269  if (priorPos == nullptr) {
270  TransientVertex tv = fitter->vertex(trTracks);
271  fittedVertex = tv;
272  } else {
273  if (priorError == nullptr) {
274  TransientVertex tv = fitter->vertex(trTracks, *priorPos);
275  fittedVertex = tv;
276  } else {
277  TransientVertex tv = fitter->vertex(trTracks, *priorPos, *priorError);
278  fittedVertex = tv;
279  }
280  }
281  } else {
282  TransientVertex tv = fitter->vertex(trTracks, *bs);
283  fittedVertex = tv;
284  }
285  } catch (std::exception const&) {
286  reco::Vertex tv;
287  fittedVertex = tv;
288  edm::LogPrint("FitFailed") << "BPHDecayVertex::fitVertex: "
289  << "vertex fit failed";
290  }
291  return;
292 }
virtual const reco::Candidate * originalReco(const reco::Candidate *daug) const
get the original particle from the clone
const BPHEventSetupWrapper * getEventSetup() const
retrieve EventSetup
const BPHEventSetupWrapper * evSetup
BPHGenericPtr< const BPHRecoCandidate >::type BPHRecoConstCandPtr
VertexFitter< 5 > * savedFitter
reco::TransientTrack * getTransientTrack(const reco::Candidate *cand) const
get TransientTrack for a daughter
virtual bool validTracks() const
check for valid reconstructed vertex
std::vector< reco::TransientTrack > trTracks
const std::vector< Component > & componentList() const
const std::vector< const reco::Track * > & tracks() const
get list of Tracks
std::map< std::string, const reco::Candidate * > dMap
std::map< const reco::Candidate *, std::string > searchMap
const reco::Track * getTrack(const reco::Candidate *cand) const
get Track for a daughter
virtual const std::vector< const reco::Candidate * > & daughFull() const
reco::TransientTrack build(const reco::Track *p) const
const std::string & getTrackSearchList(const reco::Candidate *cand) const
retrieve track search list
virtual void fitVertex(VertexFitter< 5 > *fitter, const reco::BeamSpot *bs, const GlobalPoint *priorPos, const GlobalError *priorError) const
virtual CachingVertex< N > vertex(const std::vector< reco::TransientTrack > &tracks) const =0
char getTMode(const reco::Candidate *cand) const
get Track mode for a daughter
virtual void addV(const std::string &name, const reco::Candidate *daug, const std::string &searchList, double mass)
static const reco::Track * getTrack(const reco::Candidate &rc, const char *modeList="cfhbpmnigset", char *modeFlag=nullptr)
std::map< const reco::Candidate *, reco::TransientTrack * > ttMap
virtual int charge() const =0
electric charge
Log< level::Warning, true > LogPrint
BPHDecayVertex(const BPHDecayVertex &x)=delete
virtual void addP(const std::string &name, const reco::Candidate *daug, double mass=-1.0)
reco::Vertex fittedVertex
const std::vector< reco::TransientTrack > & transientTracks() const
get list of TransientTracks
virtual const reco::Vertex & vertex(VertexFitter< 5 > *fitter=nullptr, const reco::BeamSpot *bs=nullptr, const GlobalPoint *priorPos=nullptr, const GlobalError *priorError=nullptr) const
get reconstructed vertex
virtual void fTracks() const
const GlobalError * savedPE
std::map< const reco::Candidate *, const reco::Track * > tkMap
virtual const std::vector< BPHRecoConstCandPtr > & daughComp() const
std::map< const reco::Candidate *, char > tmMap
const GlobalPoint * savedPP
virtual bool validVertex() const
virtual void setNotUpdated() const
const edm::EventSetup * get() const
void setNotUpdated() const override
virtual void fTTracks() const
virtual const std::vector< const reco::Candidate * > & daughters() const
const reco::BeamSpot * savedBS
std::vector< const reco::Track * > rTracks
~BPHDecayVertex() override
bool isValid() const
Tells whether the vertex is valid.
Definition: Vertex.h:72