CMS 3D CMS Logo

BPHRecoBuilder.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 //-------------------------------
25 
26 //---------------
27 // C++ Headers --
28 //---------------
29 using namespace std;
30 
31 //-------------------
32 // Initializations --
33 //-------------------
34 
35 //----------------
36 // Constructors --
37 //----------------
38 BPHRecoBuilder::BPHRecoBuilder(const BPHEventSetupWrapper& es) : evSetup(new BPHEventSetupWrapper(es)), minPDiff(-1.0) {
39  msList.reserve(5);
40  vsList.reserve(5);
41 }
42 
43 //--------------
44 // Destructor --
45 //--------------
47  int n = sourceList.size();
48  while (n--) {
49  delete sourceList[n]->collection;
50  delete sourceList[n];
51  }
52  int m = srCompList.size();
53  while (m--)
54  delete srCompList[m];
55  while (!compCollectList.empty()) {
56  const vector<BPHRecoConstCandPtr>* cCollection = *compCollectList.begin();
57  delete cCollection;
58  compCollectList.erase(cCollection);
59  }
60  delete evSetup;
61 }
62 
63 //--------------
64 // Operations --
65 //--------------
66 BPHRecoBuilder::BPHGenericCollection* BPHRecoBuilder::createCollection(const vector<const reco::Candidate*>& candList,
67  const string& list) {
68  return new BPHSpecificCollection<vector<const reco::Candidate*> >(candList, list);
69 }
70 
71 void BPHRecoBuilder::add(const string& name, const BPHGenericCollection* collection, double mass, double msig) {
72  BPHRecoSource* rs;
73  if (sourceId.find(name) != sourceId.end()) {
74  edm::LogPrint("TooManyParticles") << "BPHRecoBuilder::add: "
75  << "Decay product already inserted with name " << name << " , skipped";
76  return;
77  }
78  rs = new BPHRecoSource;
79  rs->name = &sourceId.insert(make_pair(name, sourceList.size())).first->first;
80  rs->collection = collection;
81  rs->selector.reserve(5);
82  rs->mass = mass;
83  rs->msig = msig;
84  sourceList.push_back(rs);
85  return;
86 }
87 
88 void BPHRecoBuilder::add(const string& name, const vector<BPHRecoConstCandPtr>& collection) {
89  BPHCompSource* cs;
90  if (srCompId.find(name) != srCompId.end()) {
91  edm::LogPrint("TooManyParticles") << "BPHRecoBuilder::add: "
92  << "Decay product already inserted with name " << name << " , skipped";
93  return;
94  }
95  cs = new BPHCompSource;
96  cs->name = &srCompId.insert(make_pair(name, srCompList.size())).first->first;
97  cs->collection = &collection;
98  srCompList.push_back(cs);
99  return;
100 }
101 
102 void BPHRecoBuilder::filter(const string& name, const BPHRecoSelect& sel) const {
103  map<string, int>::const_iterator iter = sourceId.find(name);
104  if (iter == sourceId.end())
105  return;
106  BPHRecoSource* rs = sourceList[iter->second];
107  rs->selector.push_back(&sel);
108  return;
109 }
110 
111 void BPHRecoBuilder::filter(const string& name, const BPHMomentumSelect& sel) const {
112  map<string, int>::const_iterator iter = srCompId.find(name);
113  if (iter == sourceId.end())
114  return;
115  BPHCompSource* cs = srCompList[iter->second];
116  cs->momSelector.push_back(&sel);
117  return;
118 }
119 
120 void BPHRecoBuilder::filter(const string& name, const BPHVertexSelect& sel) const {
121  map<string, int>::const_iterator iter = srCompId.find(name);
122  if (iter == sourceId.end())
123  return;
124  BPHCompSource* cs = srCompList[iter->second];
125  cs->vtxSelector.push_back(&sel);
126  return;
127 }
128 
129 void BPHRecoBuilder::filter(const string& name, const BPHFitSelect& sel) const {
130  map<string, int>::const_iterator iter = srCompId.find(name);
131  if (iter == sourceId.end())
132  return;
133  BPHCompSource* cs = srCompList[iter->second];
134  cs->fitSelector.push_back(&sel);
135  return;
136 }
137 
139  msList.push_back(&sel);
140  return;
141 }
142 
144  vsList.push_back(&sel);
145  return;
146 }
147 
149  fsList.push_back(&sel);
150  return;
151 }
152 
154  int i;
155  int n;
156  n = msList.size();
157  for (i = 0; i < n; ++i) {
158  if (!msList[i]->accept(cand))
159  return false;
160  }
161  n = vsList.size();
162  for (i = 0; i < n; ++i) {
163  if (!vsList[i]->accept(cand))
164  return false;
165  }
166  n = fsList.size();
167  for (i = 0; i < n; ++i) {
168  if (!fsList[i]->accept(cand))
169  return false;
170  }
171  return true;
172 }
173 
175  minPDiff = pMin;
176  return;
177 }
178 
179 vector<BPHRecoBuilder::ComponentSet> BPHRecoBuilder::build() const {
180  daugMap.clear();
181  compMap.clear();
182  vector<ComponentSet> candList;
183  ComponentSet compSet;
184  build(candList, compSet, sourceList.begin(), sourceList.end(), srCompList.begin(), srCompList.end());
185  return candList;
186 }
187 
189 
190 const reco::Candidate* BPHRecoBuilder::getDaug(const string& name) const {
191  map<string, const reco::Candidate*>::const_iterator iter = daugMap.find(name);
192  return (iter == daugMap.end() ? nullptr : iter->second);
193 }
194 
196  map<string, BPHRecoConstCandPtr>::const_iterator iter = compMap.find(name);
197  return (iter == compMap.end() ? nullptr : iter->second);
198 }
199 
200 bool BPHRecoBuilder::sameTrack(const reco::Candidate* lCand, const reco::Candidate* rCand, double minPDifference) {
201  const reco::Track* lrcTrack = BPHTrackReference::getFromRC(*lCand);
202  const reco::Track* rrcTrack = BPHTrackReference::getFromRC(*rCand);
203  const reco::Track* lpfTrack = BPHTrackReference::getFromPF(*lCand);
204  const reco::Track* rpfTrack = BPHTrackReference::getFromPF(*rCand);
205  if ((lrcTrack != nullptr) && ((lrcTrack == rrcTrack) || (lrcTrack == rpfTrack)))
206  return true;
207  if ((lpfTrack != nullptr) && ((lpfTrack == rrcTrack) || (lpfTrack == rpfTrack)))
208  return true;
209  reco::Candidate::Vector pDiff = (lCand->momentum() - rCand->momentum());
210  reco::Candidate::Vector pMean = (lCand->momentum() + rCand->momentum());
211  double pDMod = pDiff.mag2();
212  double pMMod = pMean.mag2();
213  if (((pDMod / pMMod) < minPDifference) && (lCand->charge() == rCand->charge()))
214  return true;
215  return false;
216 }
217 
218 void BPHRecoBuilder::build(vector<ComponentSet>& compList,
219  ComponentSet& compSet,
220  vector<BPHRecoSource*>::const_iterator r_iter,
221  vector<BPHRecoSource*>::const_iterator r_iend,
222  vector<BPHCompSource*>::const_iterator c_iter,
223  vector<BPHCompSource*>::const_iterator c_iend) const {
224  if (r_iter == r_iend) {
225  if (c_iter == c_iend) {
226  compSet.compMap = compMap;
227  compList.push_back(compSet);
228  return;
229  }
230  BPHCompSource* source = *c_iter++;
231  const vector<BPHRecoConstCandPtr>* collection = source->collection;
232  vector<const BPHMomentumSelect*> momSelector = source->momSelector;
233  vector<const BPHVertexSelect*> vtxSelector = source->vtxSelector;
234  vector<const BPHFitSelect*> fitSelector = source->fitSelector;
235  int i;
236  int j;
237  int n = collection->size();
238  int m;
239  bool skip;
240  for (i = 0; i < n; ++i) {
241  skip = false;
243  if (contained(compSet, cand))
244  continue;
245  m = momSelector.size();
246  for (j = 0; j < m; ++j) {
247  if (!momSelector[j]->accept(*cand, this)) {
248  skip = true;
249  break;
250  }
251  }
252  if (skip)
253  continue;
254  m = vtxSelector.size();
255  for (j = 0; j < m; ++j) {
256  if (!vtxSelector[j]->accept(*cand, this)) {
257  skip = true;
258  break;
259  }
260  }
261  if (skip)
262  continue;
263  m = fitSelector.size();
264  for (j = 0; j < m; ++j) {
265  if (!fitSelector[j]->accept(*cand, this)) {
266  skip = true;
267  break;
268  }
269  }
270  if (skip)
271  continue;
272  compMap[*source->name] = cand;
273  build(compList, compSet, r_iter, r_iend, c_iter, c_iend);
274  compMap.erase(*source->name);
275  }
276  return;
277  }
278  BPHRecoSource* source = *r_iter++;
279  const BPHGenericCollection* collection = source->collection;
280  vector<const BPHRecoSelect*>& selector = source->selector;
281  int i;
282  int j;
283  int n = collection->size();
284  int m = selector.size();
285  bool skip;
286  for (i = 0; i < n; ++i) {
287  const reco::Candidate& cand = collection->get(i);
288  if (contained(compSet, &cand))
289  continue;
290  skip = false;
291  for (j = 0; j < m; ++j) {
292  if (!selector[j]->accept(cand, this)) {
293  skip = true;
294  break;
295  }
296  }
297  if (skip)
298  continue;
300  comp.cand = &cand;
301  comp.mass = source->mass;
302  comp.msig = source->msig;
303  comp.searchList = collection->searchList();
304  compSet.daugMap[*source->name] = comp;
305  daugMap[*source->name] = &cand;
306  build(compList, compSet, r_iter, r_iend, c_iter, c_iend);
307  daugMap.erase(*source->name);
308  compSet.daugMap.erase(*source->name);
309  }
310  return;
311 }
312 
314  map<string, BPHDecayMomentum::Component>& dMap = compSet.daugMap;
315  map<string, BPHDecayMomentum::Component>::const_iterator d_iter;
316  map<string, BPHDecayMomentum::Component>::const_iterator d_iend = dMap.end();
317  for (d_iter = dMap.begin(); d_iter != d_iend; ++d_iter) {
318  const reco::Candidate* cChk = d_iter->second.cand;
319  if (cand == cChk)
320  return true;
321  if (sameTrack(cand, cChk))
322  return true;
323  }
324  return false;
325 }
326 
328  map<string, BPHRecoConstCandPtr>::const_iterator c_iter;
329  map<string, BPHRecoConstCandPtr>::const_iterator c_iend = compMap.end();
330  const vector<const reco::Candidate*>& dCand = cCand->daughFull();
331  int j;
332  int m = dCand.size();
333  int k;
334  int l;
335  for (j = 0; j < m; ++j) {
336  const reco::Candidate* cand = cCand->originalReco(dCand[j]);
337  map<string, BPHDecayMomentum::Component>& dMap = compSet.daugMap;
338  map<string, BPHDecayMomentum::Component>::const_iterator d_iter;
339  map<string, BPHDecayMomentum::Component>::const_iterator d_iend = dMap.end();
340  for (d_iter = dMap.begin(); d_iter != d_iend; ++d_iter) {
341  const reco::Candidate* cChk = d_iter->second.cand;
342  if (cand == cChk)
343  return true;
344  if (sameTrack(cand, cChk))
345  return true;
346  }
347 
348  for (c_iter = compMap.begin(); c_iter != c_iend; ++c_iter) {
350  BPHRecoConstCandPtr cCChk = entry.second;
351  const vector<const reco::Candidate*>& dCChk = cCChk->daughFull();
352  l = dCChk.size();
353  for (k = 0; k < l; ++k) {
354  const reco::Candidate* cChk = cCChk->originalReco(dCChk[k]);
355  if (cand == cChk)
356  return true;
357  if (sameTrack(cand, cChk))
358  return true;
359  }
360  }
361  }
362 
363  return false;
364 }
365 
366 bool BPHRecoBuilder::sameTrack(const reco::Candidate* lCand, const reco::Candidate* rCand) const {
367  return sameTrack(lCand, rCand, minPDiff);
368 }
static bool sameTrack(const reco::Candidate *lCand, const reco::Candidate *rCand, double minPDifference)
void setMinPDiffererence(double pMin)
static const reco::Track * getFromPF(const reco::Candidate &rc)
math::XYZVector Vector
point in the space
Definition: Candidate.h:42
std::map< std::string, BPHRecoConstCandPtr > compMap
bool accept(const BPHRecoCandidate &cand) const
BPHGenericPtr< const BPHRecoCandidate >::type BPHRecoConstCandPtr
virtual ~BPHRecoBuilder()
std::vector< const BPHFitSelect * > fsList
BPHRecoBuilder(const BPHEventSetupWrapper &es)
std::map< std::string, BPHDecayMomentum::Component > daugMap
const BPHEventSetupWrapper * evSetup
std::vector< BPHCompSource * > srCompList
static BPHGenericCollection * createCollection(const edm::Handle< T > &collection, const std::string &list="cfhpmig")
std::vector< const BPHVertexSelect * > vsList
static const reco::Track * getFromRC(const reco::Candidate &rc)
std::map< std::string, const reco::Candidate * > daugMap
virtual Vector momentum() const =0
spatial momentum vector
bool contained(ComponentSet &compSet, const reco::Candidate *cand) const
BPHRecoConstCandPtr getComp(const std::string &name) const
void filter(const std::string &name, const BPHRecoSelect &sel) const
virtual int charge() const =0
electric charge
std::vector< BPHRecoSource * > sourceList
Log< level::Warning, true > LogPrint
std::vector< ComponentSet > build() const
build a set of combinations of particles fulfilling the selections
void add(const std::string &name, const BPHGenericCollection *collection, double mass=-1.0, double msig=-1.0)
const reco::Candidate * getDaug(const std::string &name) const
get simple or previously recontructed particle in current combination
common object to interface with edm collections
const BPHEventSetupWrapper * eventSetup() const
get the EventSetup set in the constructor
std::map< std::string, BPHRecoConstCandPtr > compMap
std::vector< const BPHMomentumSelect * > msList
std::set< const std::vector< BPHRecoConstCandPtr > * > compCollectList
std::map< std::string, int > srCompId
static std::string const source
Definition: EdmProvDump.cc:49
std::map< std::string, int > sourceId