File indexing completed on 2024-04-06 12:15:34
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011 #include "HeavyFlavorAnalysis/RecoDecay/interface/BPHRecoBuilder.h"
0012
0013
0014
0015
0016 #include "HeavyFlavorAnalysis/RecoDecay/interface/BPHAnalyzerTokenWrapper.h"
0017 #include "HeavyFlavorAnalysis/RecoDecay/interface/BPHRecoSelect.h"
0018 #include "HeavyFlavorAnalysis/RecoDecay/interface/BPHMomentumSelect.h"
0019 #include "HeavyFlavorAnalysis/RecoDecay/interface/BPHVertexSelect.h"
0020 #include "HeavyFlavorAnalysis/RecoDecay/interface/BPHFitSelect.h"
0021 #include "HeavyFlavorAnalysis/RecoDecay/interface/BPHRecoCandidate.h"
0022 #include "HeavyFlavorAnalysis/RecoDecay/interface/BPHTrackReference.h"
0023 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
0024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0025
0026
0027
0028
0029 using namespace std;
0030
0031
0032
0033
0034
0035
0036
0037
0038 BPHRecoBuilder::BPHRecoBuilder(const BPHEventSetupWrapper& es) : evSetup(new BPHEventSetupWrapper(es)), minPDiff(-1.0) {
0039 msList.reserve(5);
0040 vsList.reserve(5);
0041 }
0042
0043
0044
0045
0046 BPHRecoBuilder::~BPHRecoBuilder() {
0047 int n = sourceList.size();
0048 while (n--) {
0049 delete sourceList[n]->collection;
0050 delete sourceList[n];
0051 }
0052 int m = srCompList.size();
0053 while (m--)
0054 delete srCompList[m];
0055 while (!compCollectList.empty()) {
0056 const vector<BPHRecoConstCandPtr>* cCollection = *compCollectList.begin();
0057 delete cCollection;
0058 compCollectList.erase(cCollection);
0059 }
0060 delete evSetup;
0061 }
0062
0063
0064
0065
0066 BPHRecoBuilder::BPHGenericCollection* BPHRecoBuilder::createCollection(const vector<const reco::Candidate*>& candList,
0067 const string& list) {
0068 return new BPHSpecificCollection<vector<const reco::Candidate*> >(candList, list);
0069 }
0070
0071 void BPHRecoBuilder::add(const string& name, const BPHGenericCollection* collection, double mass, double msig) {
0072 BPHRecoSource* rs;
0073 if (sourceId.find(name) != sourceId.end()) {
0074 edm::LogPrint("TooManyParticles") << "BPHRecoBuilder::add: "
0075 << "Decay product already inserted with name " << name << " , skipped";
0076 return;
0077 }
0078 rs = new BPHRecoSource;
0079 rs->name = &sourceId.insert(make_pair(name, sourceList.size())).first->first;
0080 rs->collection = collection;
0081 rs->selector.reserve(5);
0082 rs->mass = mass;
0083 rs->msig = msig;
0084 sourceList.push_back(rs);
0085 return;
0086 }
0087
0088 void BPHRecoBuilder::add(const string& name, const vector<BPHRecoConstCandPtr>& collection) {
0089 BPHCompSource* cs;
0090 if (srCompId.find(name) != srCompId.end()) {
0091 edm::LogPrint("TooManyParticles") << "BPHRecoBuilder::add: "
0092 << "Decay product already inserted with name " << name << " , skipped";
0093 return;
0094 }
0095 cs = new BPHCompSource;
0096 cs->name = &srCompId.insert(make_pair(name, srCompList.size())).first->first;
0097 cs->collection = &collection;
0098 srCompList.push_back(cs);
0099 return;
0100 }
0101
0102 void BPHRecoBuilder::filter(const string& name, const BPHRecoSelect& sel) const {
0103 map<string, int>::const_iterator iter = sourceId.find(name);
0104 if (iter == sourceId.end())
0105 return;
0106 BPHRecoSource* rs = sourceList[iter->second];
0107 rs->selector.push_back(&sel);
0108 return;
0109 }
0110
0111 void BPHRecoBuilder::filter(const string& name, const BPHMomentumSelect& sel) const {
0112 map<string, int>::const_iterator iter = srCompId.find(name);
0113 if (iter == sourceId.end())
0114 return;
0115 BPHCompSource* cs = srCompList[iter->second];
0116 cs->momSelector.push_back(&sel);
0117 return;
0118 }
0119
0120 void BPHRecoBuilder::filter(const string& name, const BPHVertexSelect& sel) const {
0121 map<string, int>::const_iterator iter = srCompId.find(name);
0122 if (iter == sourceId.end())
0123 return;
0124 BPHCompSource* cs = srCompList[iter->second];
0125 cs->vtxSelector.push_back(&sel);
0126 return;
0127 }
0128
0129 void BPHRecoBuilder::filter(const string& name, const BPHFitSelect& sel) const {
0130 map<string, int>::const_iterator iter = srCompId.find(name);
0131 if (iter == sourceId.end())
0132 return;
0133 BPHCompSource* cs = srCompList[iter->second];
0134 cs->fitSelector.push_back(&sel);
0135 return;
0136 }
0137
0138 void BPHRecoBuilder::filter(const BPHMomentumSelect& sel) {
0139 msList.push_back(&sel);
0140 return;
0141 }
0142
0143 void BPHRecoBuilder::filter(const BPHVertexSelect& sel) {
0144 vsList.push_back(&sel);
0145 return;
0146 }
0147
0148 void BPHRecoBuilder::filter(const BPHFitSelect& sel) {
0149 fsList.push_back(&sel);
0150 return;
0151 }
0152
0153 bool BPHRecoBuilder::accept(const BPHRecoCandidate& cand) const {
0154 int i;
0155 int n;
0156 n = msList.size();
0157 for (i = 0; i < n; ++i) {
0158 if (!msList[i]->accept(cand))
0159 return false;
0160 }
0161 n = vsList.size();
0162 for (i = 0; i < n; ++i) {
0163 if (!vsList[i]->accept(cand))
0164 return false;
0165 }
0166 n = fsList.size();
0167 for (i = 0; i < n; ++i) {
0168 if (!fsList[i]->accept(cand))
0169 return false;
0170 }
0171 return true;
0172 }
0173
0174 void BPHRecoBuilder::setMinPDiffererence(double pMin) {
0175 minPDiff = pMin;
0176 return;
0177 }
0178
0179 vector<BPHRecoBuilder::ComponentSet> BPHRecoBuilder::build() const {
0180 daugMap.clear();
0181 compMap.clear();
0182 vector<ComponentSet> candList;
0183 ComponentSet compSet;
0184 build(candList, compSet, sourceList.begin(), sourceList.end(), srCompList.begin(), srCompList.end());
0185 return candList;
0186 }
0187
0188 const BPHEventSetupWrapper* BPHRecoBuilder::eventSetup() const { return evSetup; }
0189
0190 const reco::Candidate* BPHRecoBuilder::getDaug(const string& name) const {
0191 map<string, const reco::Candidate*>::const_iterator iter = daugMap.find(name);
0192 return (iter == daugMap.end() ? nullptr : iter->second);
0193 }
0194
0195 BPHRecoConstCandPtr BPHRecoBuilder::getComp(const string& name) const {
0196 map<string, BPHRecoConstCandPtr>::const_iterator iter = compMap.find(name);
0197 return (iter == compMap.end() ? nullptr : iter->second);
0198 }
0199
0200 bool BPHRecoBuilder::sameTrack(const reco::Candidate* lCand, const reco::Candidate* rCand, double minPDifference) {
0201 const reco::Track* lrcTrack = BPHTrackReference::getFromRC(*lCand);
0202 const reco::Track* rrcTrack = BPHTrackReference::getFromRC(*rCand);
0203 const reco::Track* lpfTrack = BPHTrackReference::getFromPF(*lCand);
0204 const reco::Track* rpfTrack = BPHTrackReference::getFromPF(*rCand);
0205 if ((lrcTrack != nullptr) && ((lrcTrack == rrcTrack) || (lrcTrack == rpfTrack)))
0206 return true;
0207 if ((lpfTrack != nullptr) && ((lpfTrack == rrcTrack) || (lpfTrack == rpfTrack)))
0208 return true;
0209 reco::Candidate::Vector pDiff = (lCand->momentum() - rCand->momentum());
0210 reco::Candidate::Vector pMean = (lCand->momentum() + rCand->momentum());
0211 double pDMod = pDiff.mag2();
0212 double pMMod = pMean.mag2();
0213 if (((pDMod / pMMod) < minPDifference) && (lCand->charge() == rCand->charge()))
0214 return true;
0215 return false;
0216 }
0217
0218 void BPHRecoBuilder::build(vector<ComponentSet>& compList,
0219 ComponentSet& compSet,
0220 vector<BPHRecoSource*>::const_iterator r_iter,
0221 vector<BPHRecoSource*>::const_iterator r_iend,
0222 vector<BPHCompSource*>::const_iterator c_iter,
0223 vector<BPHCompSource*>::const_iterator c_iend) const {
0224 if (r_iter == r_iend) {
0225 if (c_iter == c_iend) {
0226 compSet.compMap = compMap;
0227 compList.push_back(compSet);
0228 return;
0229 }
0230 BPHCompSource* source = *c_iter++;
0231 const vector<BPHRecoConstCandPtr>* collection = source->collection;
0232 vector<const BPHMomentumSelect*> momSelector = source->momSelector;
0233 vector<const BPHVertexSelect*> vtxSelector = source->vtxSelector;
0234 vector<const BPHFitSelect*> fitSelector = source->fitSelector;
0235 int i;
0236 int j;
0237 int n = collection->size();
0238 int m;
0239 bool skip;
0240 for (i = 0; i < n; ++i) {
0241 skip = false;
0242 BPHRecoConstCandPtr cand = collection->at(i);
0243 if (contained(compSet, cand))
0244 continue;
0245 m = momSelector.size();
0246 for (j = 0; j < m; ++j) {
0247 if (!momSelector[j]->accept(*cand, this)) {
0248 skip = true;
0249 break;
0250 }
0251 }
0252 if (skip)
0253 continue;
0254 m = vtxSelector.size();
0255 for (j = 0; j < m; ++j) {
0256 if (!vtxSelector[j]->accept(*cand, this)) {
0257 skip = true;
0258 break;
0259 }
0260 }
0261 if (skip)
0262 continue;
0263 m = fitSelector.size();
0264 for (j = 0; j < m; ++j) {
0265 if (!fitSelector[j]->accept(*cand, this)) {
0266 skip = true;
0267 break;
0268 }
0269 }
0270 if (skip)
0271 continue;
0272 compMap[*source->name] = cand;
0273 build(compList, compSet, r_iter, r_iend, c_iter, c_iend);
0274 compMap.erase(*source->name);
0275 }
0276 return;
0277 }
0278 BPHRecoSource* source = *r_iter++;
0279 const BPHGenericCollection* collection = source->collection;
0280 vector<const BPHRecoSelect*>& selector = source->selector;
0281 int i;
0282 int j;
0283 int n = collection->size();
0284 int m = selector.size();
0285 bool skip;
0286 for (i = 0; i < n; ++i) {
0287 const reco::Candidate& cand = collection->get(i);
0288 if (contained(compSet, &cand))
0289 continue;
0290 skip = false;
0291 for (j = 0; j < m; ++j) {
0292 if (!selector[j]->accept(cand, this)) {
0293 skip = true;
0294 break;
0295 }
0296 }
0297 if (skip)
0298 continue;
0299 BPHDecayMomentum::Component comp;
0300 comp.cand = &cand;
0301 comp.mass = source->mass;
0302 comp.msig = source->msig;
0303 comp.searchList = collection->searchList();
0304 compSet.daugMap[*source->name] = comp;
0305 daugMap[*source->name] = &cand;
0306 build(compList, compSet, r_iter, r_iend, c_iter, c_iend);
0307 daugMap.erase(*source->name);
0308 compSet.daugMap.erase(*source->name);
0309 }
0310 return;
0311 }
0312
0313 bool BPHRecoBuilder::contained(ComponentSet& compSet, const reco::Candidate* cand) const {
0314 map<string, BPHDecayMomentum::Component>& dMap = compSet.daugMap;
0315 map<string, BPHDecayMomentum::Component>::const_iterator d_iter;
0316 map<string, BPHDecayMomentum::Component>::const_iterator d_iend = dMap.end();
0317 for (d_iter = dMap.begin(); d_iter != d_iend; ++d_iter) {
0318 const reco::Candidate* cChk = d_iter->second.cand;
0319 if (cand == cChk)
0320 return true;
0321 if (sameTrack(cand, cChk))
0322 return true;
0323 }
0324 return false;
0325 }
0326
0327 bool BPHRecoBuilder::contained(ComponentSet& compSet, BPHRecoConstCandPtr cCand) const {
0328 map<string, BPHRecoConstCandPtr>::const_iterator c_iter;
0329 map<string, BPHRecoConstCandPtr>::const_iterator c_iend = compMap.end();
0330 const vector<const reco::Candidate*>& dCand = cCand->daughFull();
0331 int j;
0332 int m = dCand.size();
0333 int k;
0334 int l;
0335 for (j = 0; j < m; ++j) {
0336 const reco::Candidate* cand = cCand->originalReco(dCand[j]);
0337 map<string, BPHDecayMomentum::Component>& dMap = compSet.daugMap;
0338 map<string, BPHDecayMomentum::Component>::const_iterator d_iter;
0339 map<string, BPHDecayMomentum::Component>::const_iterator d_iend = dMap.end();
0340 for (d_iter = dMap.begin(); d_iter != d_iend; ++d_iter) {
0341 const reco::Candidate* cChk = d_iter->second.cand;
0342 if (cand == cChk)
0343 return true;
0344 if (sameTrack(cand, cChk))
0345 return true;
0346 }
0347
0348 for (c_iter = compMap.begin(); c_iter != c_iend; ++c_iter) {
0349 const map<string, BPHRecoConstCandPtr>::value_type& entry = *c_iter;
0350 BPHRecoConstCandPtr cCChk = entry.second;
0351 const vector<const reco::Candidate*>& dCChk = cCChk->daughFull();
0352 l = dCChk.size();
0353 for (k = 0; k < l; ++k) {
0354 const reco::Candidate* cChk = cCChk->originalReco(dCChk[k]);
0355 if (cand == cChk)
0356 return true;
0357 if (sameTrack(cand, cChk))
0358 return true;
0359 }
0360 }
0361 }
0362
0363 return false;
0364 }
0365
0366 bool BPHRecoBuilder::sameTrack(const reco::Candidate* lCand, const reco::Candidate* rCand) const {
0367 return sameTrack(lCand, rCand, minPDiff);
0368 }