Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:05

0001 #include "DQM/Physics/src/HeavyFlavorDQMAnalyzer.h"
0002 
0003 void throwMissingCollection(const char* requested, const char* missing) {
0004   throw cms::Exception("Configuration") << "Requested plots for from the collection " << requested
0005                                         << " also requires a collection in " << missing << std::endl;
0006 }
0007 
0008 float getMass(pat::CompositeCandidate const& cand) {
0009   float mass = cand.mass();
0010   if (cand.hasUserFloat("fitMass")) {
0011     mass = cand.userFloat("fitMass");
0012   }
0013   return mass;
0014 }
0015 
0016 //
0017 // constructors and destructor
0018 //
0019 HeavyFlavorDQMAnalyzer::HeavyFlavorDQMAnalyzer(const edm::ParameterSet& iConfig)
0020     : folder_(iConfig.getParameter<std::string>("folder")),
0021       pvCollectionToken(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("pvCollection"))),
0022       beamSpotToken(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"))) {
0023   if (iConfig.existsAs<edm::InputTag>("OniaToMuMuCands")) {
0024     oniaToMuMuCandsToken =
0025         consumes<pat::CompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("OniaToMuMuCands"));
0026   }
0027   if (iConfig.existsAs<edm::InputTag>("BuToJPsiKCands")) {
0028     if (oniaToMuMuCandsToken.isUninitialized())
0029       throwMissingCollection("BuToJPsiKCands", "OniaToMuMuCands");
0030     buToJPsiKCandsToken =
0031         consumes<pat::CompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("BuToJPsiKCands"));
0032   }
0033   if (iConfig.existsAs<edm::InputTag>("Kx0ToKPiCands")) {
0034     kx0ToKPiCandsToken =
0035         consumes<pat::CompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("Kx0ToKPiCands"));
0036   }
0037   if (iConfig.existsAs<edm::InputTag>("BdToJPsiKx0Cands")) {
0038     if (oniaToMuMuCandsToken.isUninitialized())
0039       throwMissingCollection("BdToJPsiKx0Cands", "OniaToMuMuCands");
0040     if (kx0ToKPiCandsToken.isUninitialized())
0041       throwMissingCollection("BdToJPsiKx0Cands", "Kx0ToKPiCands");
0042     bdToJPsiKx0CandsToken =
0043         consumes<pat::CompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("BdToJPsiKx0Cands"));
0044   }
0045   if (iConfig.existsAs<edm::InputTag>("PhiToKKCands")) {
0046     phiToKKCandsToken =
0047         consumes<pat::CompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("PhiToKKCands"));
0048   }
0049   if (iConfig.existsAs<edm::InputTag>("BsToJPsiPhiCands")) {
0050     if (oniaToMuMuCandsToken.isUninitialized())
0051       throwMissingCollection("BsToJPsiPhiCands", "OniaToMuMuCands");
0052     if (phiToKKCandsToken.isUninitialized())
0053       throwMissingCollection("BsToJPsiPhiCands", "PhiToKKCands");
0054     bsToJPsiPhiCandsToken =
0055         consumes<pat::CompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("BsToJPsiPhiCands"));
0056   }
0057   if (iConfig.existsAs<edm::InputTag>("K0sToPiPiCands")) {
0058     k0sToPiPiCandsToken =
0059         consumes<pat::CompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("K0sToPiPiCands"));
0060   }
0061   if (iConfig.existsAs<edm::InputTag>("BdToJPsiK0sCands")) {
0062     if (oniaToMuMuCandsToken.isUninitialized())
0063       throwMissingCollection("BdToJPsiK0sCands", "OniaToMuMuCands");
0064     if (k0sToPiPiCandsToken.isUninitialized())
0065       throwMissingCollection("BdToJPsiK0sCands", "K0sToPiPiCands");
0066     bdToJPsiK0sCandsToken =
0067         consumes<pat::CompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("BdToJPsiK0sCands"));
0068   }
0069   if (iConfig.existsAs<edm::InputTag>("Lambda0ToPPiCands")) {
0070     lambda0ToPPiCandsToken =
0071         consumes<pat::CompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("Lambda0ToPPiCands"));
0072   }
0073   if (iConfig.existsAs<edm::InputTag>("LambdaBToJPsiLambda0Cands")) {
0074     if (oniaToMuMuCandsToken.isUninitialized())
0075       throwMissingCollection("LambdaBToJPsiLambda0Cands", "OniaToMuMuCands");
0076     if (lambda0ToPPiCandsToken.isUninitialized())
0077       throwMissingCollection("LambdaBToJPsiLambda0Cands", "Lambda0ToPPiCands");
0078     lambdaBToJPsiLambda0CandsToken =
0079         consumes<pat::CompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("LambdaBToJPsiLambda0Cands"));
0080   }
0081   if (iConfig.existsAs<edm::InputTag>("BcToJPsiPiCands")) {
0082     if (oniaToMuMuCandsToken.isUninitialized())
0083       throwMissingCollection("BcToJPsiPiCands", "OniaToMuMuCands");
0084     bcToJPsiPiCandsToken =
0085         consumes<pat::CompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("BcToJPsiPiCands"));
0086   }
0087   if (iConfig.existsAs<edm::InputTag>("Psi2SToJPsiPiPiCands")) {
0088     if (oniaToMuMuCandsToken.isUninitialized())
0089       throwMissingCollection("Psi2SToJPsiPiPiCands", "OniaToMuMuCands");
0090     psi2SToJPsiPiPiCandsToken =
0091         consumes<pat::CompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("Psi2SToJPsiPiPiCands"));
0092   }
0093   if (iConfig.existsAs<edm::InputTag>("BuToPsi2SKCands")) {
0094     if (psi2SToJPsiPiPiCandsToken.isUninitialized())
0095       throwMissingCollection("BuToPsi2SKCands", "Psi2SToJPsiPiPiCands");
0096     buToPsi2SKCandsToken =
0097         consumes<pat::CompositeCandidateCollection>(iConfig.getParameter<edm::InputTag>("BuToPsi2SKCands"));
0098   }
0099 }
0100 
0101 HeavyFlavorDQMAnalyzer::~HeavyFlavorDQMAnalyzer() {
0102   // do anything here that needs to be done at desctruction time
0103   // (e.g. close files, deallocate resources etc.)
0104 }
0105 
0106 //
0107 // member functions
0108 //
0109 
0110 // ------------ method called for each event  ------------
0111 
0112 void HeavyFlavorDQMAnalyzer::dqmAnalyze(edm::Event const& iEvent,
0113                                         edm::EventSetup const& iSetup,
0114                                         Histograms const& histos) const {
0115   auto& pvColl = iEvent.get(pvCollectionToken);
0116   auto bs = iEvent.getHandle(beamSpotToken).product();
0117 
0118   std::vector<bool> displacedJPsiToMuMu;
0119   pat::CompositeCandidateCollection lOniaToMuMu;
0120   if (not oniaToMuMuCandsToken.isUninitialized()) {
0121     lOniaToMuMu = iEvent.get(oniaToMuMuCandsToken);
0122     displacedJPsiToMuMu.resize(lOniaToMuMu.size(), false);
0123   }
0124 
0125   if (not buToJPsiKCandsToken.isUninitialized()) {
0126     auto lBuToJPsiK = iEvent.get(buToJPsiKCandsToken);
0127     for (auto&& cand : lBuToJPsiK) {
0128       auto jpsi = *cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi");
0129       auto jpsiMass = getMass(*jpsi);
0130       if (jpsiMass < 2.9 or jpsiMass > 3.3)
0131         continue;
0132 
0133       auto closestPV = fillDecayHistograms(histos.buToJPsiK, cand, pvColl);
0134       if (not closestPV)
0135         continue;
0136       fillBuToJPsiKComponents(histos.buToJPsiK, cand, bs, closestPV);
0137 
0138       displacedJPsiToMuMu[jpsi.index()] = true;
0139     }
0140   }
0141 
0142   std::vector<bool> displacedPsi2SToJPsiPiPi;
0143   pat::CompositeCandidateCollection lPsi2SToJPsiPiPi;
0144   if (not psi2SToJPsiPiPiCandsToken.isUninitialized()) {
0145     lPsi2SToJPsiPiPi = iEvent.get(psi2SToJPsiPiPiCandsToken);
0146     displacedPsi2SToJPsiPiPi.resize(lPsi2SToJPsiPiPi.size(), false);
0147   }
0148 
0149   if (not buToPsi2SKCandsToken.isUninitialized()) {
0150     auto lBuToPsi2SK = iEvent.get(buToPsi2SKCandsToken);
0151     for (auto&& cand : lBuToPsi2SK) {
0152       auto psi2S = *cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToPsi2S");
0153       auto psi2SMass = getMass(*psi2S);
0154       if (psi2SMass < 3.65 or psi2SMass > 3.72)
0155         continue;
0156 
0157       auto jpsi = *psi2S->userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi");
0158       auto jpsiMass = getMass(*jpsi);
0159       if (jpsiMass < 2.9 or jpsiMass > 3.3)
0160         continue;
0161 
0162       auto closestPV = fillDecayHistograms(histos.buToPsi2SK, cand, pvColl);
0163       if (not closestPV)
0164         continue;
0165       fillBuToPsi2SKComponents(histos.buToPsi2SK, cand, bs, closestPV);
0166 
0167       displacedPsi2SToJPsiPiPi[psi2S.index()] = true;
0168       displacedJPsiToMuMu[jpsi.index()] = true;
0169     }
0170   }
0171 
0172   for (size_t i = 0; i < lPsi2SToJPsiPiPi.size(); i++) {
0173     auto&& cand = lPsi2SToJPsiPiPi[i];
0174 
0175     auto jpsi = *cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi");
0176     auto jpsiMass = getMass(*jpsi);
0177     if (jpsiMass < 2.9 or jpsiMass > 3.3)
0178       continue;
0179 
0180     auto closestPV = fillDecayHistograms(histos.psi2SToJPsiPiPi, cand, pvColl);
0181     if (not closestPV)
0182       continue;
0183     fillPsi2SToJPsiPiPiComponents(histos.psi2SToJPsiPiPi, cand, bs, closestPV);
0184 
0185     auto decayHistos = &histos.psi2SToJPsiPiPiPrompt;
0186     if (displacedPsi2SToJPsiPiPi[i]) {
0187       decayHistos = &histos.psi2SToJPsiPiPiDispl;
0188     }
0189 
0190     fillDecayHistograms(*decayHistos, cand, pvColl);
0191     fillPsi2SToJPsiPiPiComponents(*decayHistos, cand, bs, closestPV);
0192   }
0193   lPsi2SToJPsiPiPi.clear();
0194   displacedPsi2SToJPsiPiPi.clear();
0195 
0196   std::vector<bool> displacedKx0ToKPi;
0197   pat::CompositeCandidateCollection lKx0ToKPi;
0198   if (not kx0ToKPiCandsToken.isUninitialized()) {
0199     lKx0ToKPi = iEvent.get(kx0ToKPiCandsToken);
0200     displacedKx0ToKPi.resize(lKx0ToKPi.size(), false);
0201   }
0202 
0203   if (not bdToJPsiKx0CandsToken.isUninitialized()) {
0204     auto lBdToJPsiKx0 = iEvent.get(bdToJPsiKx0CandsToken);
0205     for (auto&& cand : lBdToJPsiKx0) {
0206       auto jpsi = *cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi");
0207       auto jpsiMass = getMass(*jpsi);
0208       if (jpsiMass < 2.9 or jpsiMass > 3.3)
0209         continue;
0210 
0211       auto kx0 = *cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToKx0");
0212       auto kx0Mass = getMass(*kx0);
0213       if (kx0Mass < 0.77 or kx0Mass > 1.02)
0214         continue;
0215 
0216       auto closestPV = fillDecayHistograms(histos.bdToJPsiKx0, cand, pvColl);
0217       if (not closestPV)
0218         continue;
0219       fillBdToJPsiKx0Components(histos.bdToJPsiKx0, cand, bs, closestPV);
0220 
0221       displacedKx0ToKPi[kx0.index()] = true;
0222     }
0223   }
0224 
0225   for (size_t i = 0; i < lKx0ToKPi.size(); i++) {
0226     auto&& cand = lKx0ToKPi[i];
0227 
0228     auto closestPV = fillDecayHistograms(histos.kx0ToKPi, cand, pvColl);
0229     if (not closestPV)
0230       continue;
0231     fillKx0ToKPiComponents(histos.kx0ToKPi, cand, bs, closestPV);
0232 
0233     auto decayHistos = &histos.kx0ToKPiPrompt;
0234     if (displacedKx0ToKPi[i]) {
0235       decayHistos = &histos.kx0ToKPiDispl;
0236     }
0237 
0238     fillDecayHistograms(*decayHistos, cand, pvColl);
0239     fillKx0ToKPiComponents(*decayHistos, cand, bs, closestPV);
0240   }
0241   lKx0ToKPi.clear();
0242   displacedKx0ToKPi.clear();
0243 
0244   std::vector<bool> displacedPhiToKK;
0245   pat::CompositeCandidateCollection lPhiToKK;
0246   if (not phiToKKCandsToken.isUninitialized()) {
0247     lPhiToKK = iEvent.get(phiToKKCandsToken);
0248     displacedPhiToKK.resize(lPhiToKK.size(), false);
0249   }
0250 
0251   if (not bsToJPsiPhiCandsToken.isUninitialized()) {
0252     auto lBsToJPsiPhi = iEvent.get(bsToJPsiPhiCandsToken);
0253     for (auto&& cand : lBsToJPsiPhi) {
0254       auto jpsi = *cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi");
0255       auto jpsiMass = getMass(*jpsi);
0256       if (jpsiMass < 2.9 or jpsiMass > 3.3)
0257         continue;
0258 
0259       auto phi = *cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToPhi");
0260       auto phiMass = getMass(*phi);
0261       if (phiMass < 1.005 or phiMass > 1.035)
0262         continue;
0263 
0264       auto closestPV = fillDecayHistograms(histos.bsToJPsiPhi, cand, pvColl);
0265       if (not closestPV)
0266         continue;
0267       fillBsToJPsiPhiComponents(histos.bsToJPsiPhi, cand, bs, closestPV);
0268 
0269       displacedJPsiToMuMu[jpsi.index()] = true;
0270       displacedPhiToKK[phi.index()] = true;
0271     }
0272   }
0273 
0274   for (size_t i = 0; i < lPhiToKK.size(); i++) {
0275     auto&& cand = lPhiToKK[i];
0276 
0277     auto closestPV = fillDecayHistograms(histos.phiToKK, cand, pvColl);
0278     if (not closestPV)
0279       continue;
0280     fillPhiToKKComponents(histos.phiToKK, cand, bs, closestPV);
0281 
0282     auto decayHistos = &histos.phiToKKPrompt;
0283     if (displacedPhiToKK[i]) {
0284       decayHistos = &histos.phiToKKDispl;
0285     }
0286 
0287     fillDecayHistograms(*decayHistos, cand, pvColl);
0288     fillPhiToKKComponents(*decayHistos, cand, bs, closestPV);
0289   }
0290   lPhiToKK.clear();
0291   displacedPhiToKK.clear();
0292 
0293   if (not bdToJPsiK0sCandsToken.isUninitialized()) {
0294     auto lBdToJPsiK0s = iEvent.get(bdToJPsiK0sCandsToken);
0295     for (auto&& cand : lBdToJPsiK0s) {
0296       auto jpsi = *cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi");
0297       auto jpsiMass = getMass(*jpsi);
0298       if (jpsiMass < 2.9 or jpsiMass > 3.3)
0299         continue;
0300 
0301       auto closestPV = fillDecayHistograms(histos.bdToJPsiK0s, cand, pvColl);
0302       if (not closestPV)
0303         continue;
0304       fillBdToJPsiK0sComponents(histos.bdToJPsiK0s, cand, bs, closestPV);
0305 
0306       displacedJPsiToMuMu[jpsi.index()] = true;
0307     }
0308   }
0309 
0310   if (not bcToJPsiPiCandsToken.isUninitialized()) {
0311     auto lBcToJPsiPi = iEvent.get(bcToJPsiPiCandsToken);
0312     for (auto&& cand : lBcToJPsiPi) {
0313       auto jpsi = *cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi");
0314       auto jpsiMass = getMass(*jpsi);
0315       if (jpsiMass < 2.9 or jpsiMass > 3.3)
0316         continue;
0317 
0318       auto closestPV = fillDecayHistograms(histos.bcToJPsiPi, cand, pvColl);
0319       if (not closestPV)
0320         continue;
0321       fillBcToJPsiPiComponents(histos.bcToJPsiPi, cand, bs, closestPV);
0322 
0323       displacedJPsiToMuMu[jpsi.index()] = true;
0324     }
0325   }
0326 
0327   if (not lambdaBToJPsiLambda0CandsToken.isUninitialized()) {
0328     auto lLambdaBToJPsiLambda0 = iEvent.get(lambdaBToJPsiLambda0CandsToken);
0329     for (auto&& cand : lLambdaBToJPsiLambda0) {
0330       auto jpsi = *cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi");
0331       auto jpsiMass = getMass(*jpsi);
0332       if (jpsiMass < 2.9 or jpsiMass > 3.3)
0333         continue;
0334 
0335       auto closestPV = fillDecayHistograms(histos.lambdaBToJPsiLambda0, cand, pvColl);
0336       if (not closestPV)
0337         continue;
0338       fillLambdaBToJPsiLambda0Components(histos.lambdaBToJPsiLambda0, cand, bs, closestPV);
0339 
0340       displacedJPsiToMuMu[jpsi.index()] = true;
0341     }
0342   }
0343 
0344   for (size_t i = 0; i < lOniaToMuMu.size(); i++) {
0345     auto&& cand = lOniaToMuMu[i];
0346 
0347     auto closestPV = fillDecayHistograms(histos.oniaToMuMu, cand, pvColl);
0348     if (not closestPV)
0349       continue;
0350     fillOniaToMuMuComponents(histos.oniaToMuMu, cand, bs, closestPV);
0351 
0352     auto decayHistos = &histos.oniaToMuMuPrompt;
0353     if (displacedJPsiToMuMu[i]) {
0354       decayHistos = &histos.oniaToMuMuDispl;
0355     }
0356 
0357     fillDecayHistograms(*decayHistos, cand, pvColl);
0358     fillOniaToMuMuComponents(*decayHistos, cand, bs, closestPV);
0359   }
0360   lOniaToMuMu.clear();
0361   displacedJPsiToMuMu.clear();
0362 
0363   if (not k0sToPiPiCandsToken.isUninitialized()) {
0364     auto lK0sToPiPi = iEvent.get(k0sToPiPiCandsToken);
0365     for (auto&& cand : lK0sToPiPi) {
0366       auto closestPV = fillDecayHistograms(histos.k0sToPiPi, cand, pvColl);
0367       if (not closestPV)
0368         continue;
0369       fillK0sToPiPiComponents(histos.k0sToPiPi, cand, bs, closestPV);
0370     }
0371   }
0372 
0373   if (not lambda0ToPPiCandsToken.isUninitialized()) {
0374     auto lLambda0ToPPi = iEvent.get(lambda0ToPPiCandsToken);
0375     for (auto&& cand : lLambda0ToPPi) {
0376       auto closestPV = fillDecayHistograms(histos.lambda0ToPPi, cand, pvColl);
0377       if (not closestPV)
0378         continue;
0379       fillLambda0ToPPiComponents(histos.lambda0ToPPi, cand, bs, closestPV);
0380     }
0381   }
0382 }
0383 
0384 void HeavyFlavorDQMAnalyzer::bookDecayHists(DQMStore::IBooker& ibook,
0385                                             edm::Run const&,
0386                                             edm::EventSetup const&,
0387                                             DecayHists& decayHists,
0388                                             std::string const& name,
0389                                             std::string const& products,
0390                                             int nMassBins,
0391                                             float massMin,
0392                                             float massMax,
0393                                             float distanceScaleFactor) const {
0394   std::string histTitle = name + " #rightarrow " + products + ";";
0395 
0396   decayHists.h_mass =
0397       ibook.book1D("h_mass", histTitle + "M(" + products + ") fitted [GeV]", nMassBins, massMin, massMax);
0398   decayHists.h_pt = ibook.book1D("h_pt", histTitle + "fitted p_{T} [GeV]", 100, 0.00, 200.0);
0399   decayHists.h_eta = ibook.book1D("h_eta", histTitle + "fitted #eta", 100, -3, 3);
0400   decayHists.h_phi = ibook.book1D("h_phi", histTitle + "fitted #varphi [rad]", 100, -TMath::Pi(), TMath::Pi());
0401   decayHists.h_displ2D =
0402       ibook.book1D("h_displ2D", histTitle + "vertex 2D displacement [cm]", 100, 0.00, 2.0 * distanceScaleFactor);
0403   decayHists.h_sign2D =
0404       ibook.book1D("h_sign2D", histTitle + "vertex 2D displ. significance", 100, 0.00, 200.0 * distanceScaleFactor);
0405   decayHists.h_ct = ibook.book1D("h_ct", histTitle + "ct [cm]", 100, 0.00, 0.4 * distanceScaleFactor);
0406   decayHists.h_pointing = ibook.book1D("h_pointing", histTitle + "cos( 2D pointing angle )", 100, -1, 1);
0407   decayHists.h_vertNormChi2 = ibook.book1D("h_vertNormChi2", histTitle + "vertex #chi^{2}/ndof", 100, 0.00, 10);
0408   decayHists.h_vertProb = ibook.book1D("h_vertProb", histTitle + "vertex prob.", 100, 0.00, 1.0);
0409 }
0410 
0411 void HeavyFlavorDQMAnalyzer::bookHistograms(DQMStore::IBooker& ibook,
0412                                             edm::Run const& run,
0413                                             edm::EventSetup const& iSetup,
0414                                             Histograms& histos) const {
0415   if (not oniaToMuMuCandsToken.isUninitialized()) {
0416     ibook.cd();
0417     ibook.setCurrentFolder(folder_ + "/JPsiToMuMuPrompt");
0418     bookDecayHists(ibook, run, iSetup, histos.oniaToMuMuPrompt, "J/#psi", "#mu^{+}#mu^{-}", 100, 2.9, 3.3);
0419 
0420     ibook.setCurrentFolder(folder_ + "/JPsiToMuMuPrompt/components");
0421     initOniaToMuMuComponentHistograms(ibook, run, iSetup, histos.oniaToMuMuPrompt);
0422 
0423     ibook.cd();
0424     ibook.setCurrentFolder(folder_ + "/JPsiToMuMuDisplaced");
0425     bookDecayHists(ibook, run, iSetup, histos.oniaToMuMuDispl, "J/#psi", "#mu^{+}#mu^{-}", 100, 2.9, 3.3);
0426 
0427     ibook.setCurrentFolder(folder_ + "/JPsiToMuMuDisplaced/components");
0428     initOniaToMuMuComponentHistograms(ibook, run, iSetup, histos.oniaToMuMuDispl);
0429 
0430     ibook.cd();
0431     ibook.setCurrentFolder(folder_ + "/OniaToMuMu");
0432     bookDecayHists(ibook, run, iSetup, histos.oniaToMuMu, "Onia", "#mu^{+}#mu^{-}", 750, 0, 15);
0433 
0434     ibook.setCurrentFolder(folder_ + "/OniaToMuMu/components");
0435     initOniaToMuMuComponentHistograms(ibook, run, iSetup, histos.oniaToMuMu);
0436   }
0437 
0438   if (not kx0ToKPiCandsToken.isUninitialized()) {
0439     ibook.cd();
0440     ibook.setCurrentFolder(folder_ + "/Kx0ToKPiPrompt");
0441     bookDecayHists(ibook, run, iSetup, histos.kx0ToKPiPrompt, "K*^{0}", "#pi^{+} K^{-}", 100, 0.75, 1.05);
0442 
0443     ibook.setCurrentFolder(folder_ + "/Kx0ToKPiPrompt/components");
0444     initKx0ToKPiComponentHistograms(ibook, run, iSetup, histos.kx0ToKPiPrompt);
0445 
0446     ibook.cd();
0447     ibook.setCurrentFolder(folder_ + "/Kx0ToKPiDisplaced");
0448     bookDecayHists(ibook, run, iSetup, histos.kx0ToKPiDispl, "K*^{0}", "#pi^{+} K^{-}", 100, 0.75, 1.05);
0449 
0450     ibook.setCurrentFolder(folder_ + "/Kx0ToKPiDisplaced/components");
0451     initKx0ToKPiComponentHistograms(ibook, run, iSetup, histos.kx0ToKPiDispl);
0452 
0453     ibook.cd();
0454     ibook.setCurrentFolder(folder_ + "/Kx0ToKPi");
0455     bookDecayHists(ibook, run, iSetup, histos.kx0ToKPi, "K*^{0}", "#pi^{+} K^{-}", 100, 0.75, 1.05);
0456 
0457     ibook.setCurrentFolder(folder_ + "/Kx0ToKPi/components");
0458     initKx0ToKPiComponentHistograms(ibook, run, iSetup, histos.kx0ToKPi);
0459   }
0460 
0461   if (not phiToKKCandsToken.isUninitialized()) {
0462     ibook.cd();
0463     ibook.setCurrentFolder(folder_ + "/PhiToKKPrompt");
0464     bookDecayHists(ibook, run, iSetup, histos.phiToKKPrompt, "#phi", "K^{+} K^{-}", 100, 1.005, 1.035);
0465 
0466     ibook.setCurrentFolder(folder_ + "/PhiToKKPrompt/components");
0467     initPhiToKKComponentHistograms(ibook, run, iSetup, histos.phiToKKPrompt);
0468 
0469     ibook.cd();
0470     ibook.setCurrentFolder(folder_ + "/PhiToKKDisplaced");
0471     bookDecayHists(ibook, run, iSetup, histos.phiToKKDispl, "#phi", "K^{+} K^{-}", 100, 1.005, 1.035);
0472 
0473     ibook.setCurrentFolder(folder_ + "/PhiToKKDisplaced/components");
0474     initPhiToKKComponentHistograms(ibook, run, iSetup, histos.phiToKKDispl);
0475 
0476     ibook.cd();
0477     ibook.setCurrentFolder(folder_ + "/PhiToKK");
0478     bookDecayHists(ibook, run, iSetup, histos.phiToKK, "#phi", "K^{+} K^{-}", 100, 1.005, 1.035);
0479 
0480     ibook.setCurrentFolder(folder_ + "/PhiToKK/components");
0481     initPhiToKKComponentHistograms(ibook, run, iSetup, histos.phiToKK);
0482   }
0483 
0484   if (not psi2SToJPsiPiPiCandsToken.isUninitialized()) {
0485     ibook.cd();
0486     ibook.setCurrentFolder(folder_ + "/Psi2SToJPsiPiPiPrompt");
0487     bookDecayHists(
0488         ibook, run, iSetup, histos.psi2SToJPsiPiPiPrompt, "#Psi(2S)", "J/#psi #pi^{+} #pi^{-}", 100, 3.65, 3.72);
0489 
0490     ibook.setCurrentFolder(folder_ + "/Psi2SToJPsiPiPiPrompt/components");
0491     initPsi2SToJPsiPiPiComponentHistograms(ibook, run, iSetup, histos.psi2SToJPsiPiPiPrompt);
0492 
0493     ibook.cd();
0494     ibook.setCurrentFolder(folder_ + "/Psi2SToJPsiPiPiDisplaced");
0495     bookDecayHists(
0496         ibook, run, iSetup, histos.psi2SToJPsiPiPiDispl, "#Psi(2S)", "J/#psi #pi^{+} #pi^{-}", 100, 3.65, 3.72);
0497 
0498     ibook.setCurrentFolder(folder_ + "/Psi2SToJPsiPiPiDisplaced/components");
0499     initPsi2SToJPsiPiPiComponentHistograms(ibook, run, iSetup, histos.psi2SToJPsiPiPiDispl);
0500 
0501     ibook.cd();
0502     ibook.setCurrentFolder(folder_ + "/Psi2SToJPsiPiPi");
0503     bookDecayHists(
0504         ibook, run, iSetup, histos.psi2SToJPsiPiPi, "#Psi(2S)/X(3872)", "J/#psi #pi^{+} #pi^{-}", 200, 3.60, 3.80);
0505 
0506     ibook.setCurrentFolder(folder_ + "/Psi2SToJPsiPiPi/components");
0507     initPsi2SToJPsiPiPiComponentHistograms(ibook, run, iSetup, histos.psi2SToJPsiPiPi);
0508   }
0509 
0510   if (not k0sToPiPiCandsToken.isUninitialized()) {
0511     ibook.cd();
0512     ibook.setCurrentFolder(folder_ + "/K0sToPiPi");
0513     bookDecayHists(ibook, run, iSetup, histos.k0sToPiPi, "K^{0}_{S}", "#pi^{+} #pi^{-}", 100, 0.44, 0.56, 4);
0514 
0515     ibook.setCurrentFolder(folder_ + "/K0sToPiPi/components");
0516     initK0sToPiPiComponentHistograms(ibook, run, iSetup, histos.k0sToPiPi);
0517   }
0518 
0519   if (not lambda0ToPPiCandsToken.isUninitialized()) {
0520     ibook.cd();
0521     ibook.setCurrentFolder(folder_ + "/Lambda0ToPPi");
0522     bookDecayHists(ibook, run, iSetup, histos.lambda0ToPPi, "#Lambda^{0}", "p^{+} #pi^{-}", 100, 1.06, 1.16, 4);
0523 
0524     ibook.setCurrentFolder(folder_ + "/Lambda0ToPPi/components");
0525     initLambda0ToPPiComponentHistograms(ibook, run, iSetup, histos.lambda0ToPPi);
0526   }
0527 
0528   if (not buToJPsiKCandsToken.isUninitialized()) {
0529     ibook.cd();
0530     ibook.setCurrentFolder(folder_ + "/BuToJPsiK");
0531     bookDecayHists(ibook, run, iSetup, histos.buToJPsiK, "B^{+}", "J/#psi K^{+}", 100, 5.00, 6.00);
0532 
0533     ibook.setCurrentFolder(folder_ + "/BuToJPsiK/components");
0534     initBuToJPsiKComponentHistograms(ibook, run, iSetup, histos.buToJPsiK);
0535   }
0536 
0537   if (not buToPsi2SKCandsToken.isUninitialized()) {
0538     ibook.cd();
0539     ibook.setCurrentFolder(folder_ + "/BuToPsi2SK");
0540     bookDecayHists(ibook, run, iSetup, histos.buToPsi2SK, "B^{+}", "#Psi(2S) K^{+}", 100, 5.00, 6.00);
0541 
0542     ibook.setCurrentFolder(folder_ + "/BuToPsi2SK/components");
0543     initBuToPsi2SKComponentHistograms(ibook, run, iSetup, histos.buToPsi2SK);
0544   }
0545 
0546   if (not bdToJPsiKx0CandsToken.isUninitialized()) {
0547     ibook.cd();
0548     ibook.setCurrentFolder(folder_ + "/BdToJPsiKx0");
0549     bookDecayHists(ibook, run, iSetup, histos.bdToJPsiKx0, "B^{0}", "J/#psi K*^{0}", 100, 5.00, 6.00);
0550 
0551     ibook.setCurrentFolder(folder_ + "/BdToJPsiKx0/components");
0552     initBdToJPsiKx0ComponentHistograms(ibook, run, iSetup, histos.bdToJPsiKx0);
0553   }
0554 
0555   if (not bsToJPsiPhiCandsToken.isUninitialized()) {
0556     ibook.cd();
0557     ibook.setCurrentFolder(folder_ + "/BsToJPsiPhi");
0558     bookDecayHists(ibook, run, iSetup, histos.bsToJPsiPhi, "B^{0}_{s}", "J/#psi #phi", 100, 5.00, 6.00);
0559 
0560     ibook.setCurrentFolder(folder_ + "/BsToJPsiPhi/components");
0561     initBsToJPsiPhiComponentHistograms(ibook, run, iSetup, histos.bsToJPsiPhi);
0562   }
0563 
0564   if (not bdToJPsiK0sCandsToken.isUninitialized()) {
0565     ibook.cd();
0566     ibook.setCurrentFolder(folder_ + "/BdToJPsiK0s");
0567     bookDecayHists(ibook, run, iSetup, histos.bdToJPsiK0s, "B^{0}", "J/#psi K^{0}_{S}", 100, 5.00, 6.00);
0568 
0569     ibook.setCurrentFolder(folder_ + "/BdToJPsiK0s/components");
0570     initBdToJPsiK0sComponentHistograms(ibook, run, iSetup, histos.bdToJPsiK0s);
0571   }
0572 
0573   if (not bcToJPsiPiCandsToken.isUninitialized()) {
0574     ibook.cd();
0575     ibook.setCurrentFolder(folder_ + "/BcToJPsiPi");
0576     bookDecayHists(ibook, run, iSetup, histos.bcToJPsiPi, "B^{+}_{c}", "J/#psi #pi^{+}", 100, 6.00, 7.00);
0577 
0578     ibook.setCurrentFolder(folder_ + "/BcToJPsiPi/components");
0579     initBcToJPsiPiComponentHistograms(ibook, run, iSetup, histos.bcToJPsiPi);
0580   }
0581 
0582   if (not lambdaBToJPsiLambda0CandsToken.isUninitialized()) {
0583     ibook.cd();
0584     ibook.setCurrentFolder(folder_ + "/LambdaBToJPsiLambda0");
0585     bookDecayHists(
0586         ibook, run, iSetup, histos.lambdaBToJPsiLambda0, "#Lambda^{0}_{b}", "J/#psi #Lambda^{0}", 100, 5.00, 6.00);
0587 
0588     ibook.setCurrentFolder(folder_ + "/LambdaBToJPsiLambda0/components");
0589     initLambdaBToJPsiLambda0ComponentHistograms(ibook, run, iSetup, histos.lambdaBToJPsiLambda0);
0590   }
0591 }
0592 
0593 void HeavyFlavorDQMAnalyzer::initComponentHists(DQMStore::IBooker& ibook,
0594                                                 edm::Run const&,
0595                                                 edm::EventSetup const&,
0596                                                 DecayHists& histos,
0597                                                 TString const& componentName) const {
0598   ComponentHists comp;
0599 
0600   comp.h_pt = ibook.book1D(componentName + "_pt", ";p_{T} [GeV]", 200, 0, 20);
0601   comp.h_eta = ibook.book1D(componentName + "_eta", ";#eta", 200, -3, 3);
0602   comp.h_phi = ibook.book1D(componentName + "_phi", ";#phi", 200, -TMath::Pi(), TMath::Pi());
0603   comp.h_dxy = ibook.book1D(componentName + "_dxyBS", ";d_{xy}(BS) [cm]", 200, -3, 3);
0604   comp.h_exy = ibook.book1D(componentName + "_exy", ";#sigma(d_{xy}(BS)) [cm]", 200, 0, 0.2);
0605   comp.h_dz = ibook.book1D(componentName + "_dzPV", ";d_{z}(PV) [cm]", 200, -20, 20);
0606   comp.h_ez = ibook.book1D(componentName + "_ez", ";#sigma(d_{z}(PV)) [cm]", 200, 0, 2);
0607   comp.h_chi2 = ibook.book1D(componentName + "_chi2", ";#chi^{2}", 200, 0, 20);
0608 
0609   histos.decayComponents.push_back(comp);
0610 }
0611 
0612 void HeavyFlavorDQMAnalyzer::initOniaToMuMuComponentHistograms(DQMStore::IBooker& ibook,
0613                                                                edm::Run const& run,
0614                                                                edm::EventSetup const& iSetup,
0615                                                                DecayHists& histos) const {
0616   initComponentHists(ibook, run, iSetup, histos, "lead_mu");
0617   initComponentHists(ibook, run, iSetup, histos, "soft_mu");
0618 }
0619 
0620 void HeavyFlavorDQMAnalyzer::initKx0ToKPiComponentHistograms(DQMStore::IBooker& ibook,
0621                                                              edm::Run const& run,
0622                                                              edm::EventSetup const& iSetup,
0623                                                              DecayHists& histos) const {
0624   initComponentHists(ibook, run, iSetup, histos, "k");
0625   initComponentHists(ibook, run, iSetup, histos, "pi");
0626 }
0627 
0628 void HeavyFlavorDQMAnalyzer::initPhiToKKComponentHistograms(DQMStore::IBooker& ibook,
0629                                                             edm::Run const& run,
0630                                                             edm::EventSetup const& iSetup,
0631                                                             DecayHists& histos) const {
0632   initComponentHists(ibook, run, iSetup, histos, "lead_k");
0633   initComponentHists(ibook, run, iSetup, histos, "soft_k");
0634 }
0635 
0636 void HeavyFlavorDQMAnalyzer::initK0sToPiPiComponentHistograms(DQMStore::IBooker& ibook,
0637                                                               edm::Run const& run,
0638                                                               edm::EventSetup const& iSetup,
0639                                                               DecayHists& histos) const {
0640   initComponentHists(ibook, run, iSetup, histos, "lead_pi");
0641   initComponentHists(ibook, run, iSetup, histos, "soft_pi");
0642 }
0643 
0644 void HeavyFlavorDQMAnalyzer::initLambda0ToPPiComponentHistograms(DQMStore::IBooker& ibook,
0645                                                                  edm::Run const& run,
0646                                                                  edm::EventSetup const& iSetup,
0647                                                                  DecayHists& histos) const {
0648   initComponentHists(ibook, run, iSetup, histos, "p");
0649   initComponentHists(ibook, run, iSetup, histos, "pi");
0650 }
0651 
0652 void HeavyFlavorDQMAnalyzer::initBuToJPsiKComponentHistograms(DQMStore::IBooker& ibook,
0653                                                               edm::Run const& run,
0654                                                               edm::EventSetup const& iSetup,
0655                                                               DecayHists& histos) const {
0656   initOniaToMuMuComponentHistograms(ibook, run, iSetup, histos);
0657   initComponentHists(ibook, run, iSetup, histos, "k");
0658 }
0659 
0660 void HeavyFlavorDQMAnalyzer::initBuToPsi2SKComponentHistograms(DQMStore::IBooker& ibook,
0661                                                                edm::Run const& run,
0662                                                                edm::EventSetup const& iSetup,
0663                                                                DecayHists& histos) const {
0664   initPsi2SToJPsiPiPiComponentHistograms(ibook, run, iSetup, histos);
0665   initComponentHists(ibook, run, iSetup, histos, "k");
0666 }
0667 
0668 void HeavyFlavorDQMAnalyzer::initBdToJPsiKx0ComponentHistograms(DQMStore::IBooker& ibook,
0669                                                                 edm::Run const& run,
0670                                                                 edm::EventSetup const& iSetup,
0671                                                                 DecayHists& histos) const {
0672   initOniaToMuMuComponentHistograms(ibook, run, iSetup, histos);
0673   initKx0ToKPiComponentHistograms(ibook, run, iSetup, histos);
0674 }
0675 
0676 void HeavyFlavorDQMAnalyzer::initBsToJPsiPhiComponentHistograms(DQMStore::IBooker& ibook,
0677                                                                 edm::Run const& run,
0678                                                                 edm::EventSetup const& iSetup,
0679                                                                 DecayHists& histos) const {
0680   initOniaToMuMuComponentHistograms(ibook, run, iSetup, histos);
0681   initPhiToKKComponentHistograms(ibook, run, iSetup, histos);
0682 }
0683 
0684 void HeavyFlavorDQMAnalyzer::initBdToJPsiK0sComponentHistograms(DQMStore::IBooker& ibook,
0685                                                                 edm::Run const& run,
0686                                                                 edm::EventSetup const& iSetup,
0687                                                                 DecayHists& histos) const {
0688   initOniaToMuMuComponentHistograms(ibook, run, iSetup, histos);
0689   initK0sToPiPiComponentHistograms(ibook, run, iSetup, histos);
0690 }
0691 
0692 void HeavyFlavorDQMAnalyzer::initBcToJPsiPiComponentHistograms(DQMStore::IBooker& ibook,
0693                                                                edm::Run const& run,
0694                                                                edm::EventSetup const& iSetup,
0695                                                                DecayHists& histos) const {
0696   initOniaToMuMuComponentHistograms(ibook, run, iSetup, histos);
0697   initComponentHists(ibook, run, iSetup, histos, "pi");
0698 }
0699 
0700 void HeavyFlavorDQMAnalyzer::initLambdaBToJPsiLambda0ComponentHistograms(DQMStore::IBooker& ibook,
0701                                                                          edm::Run const& run,
0702                                                                          edm::EventSetup const& iSetup,
0703                                                                          DecayHists& histos) const {
0704   initOniaToMuMuComponentHistograms(ibook, run, iSetup, histos);
0705   initLambda0ToPPiComponentHistograms(ibook, run, iSetup, histos);
0706 }
0707 
0708 void HeavyFlavorDQMAnalyzer::initPsi2SToJPsiPiPiComponentHistograms(DQMStore::IBooker& ibook,
0709                                                                     edm::Run const& run,
0710                                                                     edm::EventSetup const& iSetup,
0711                                                                     DecayHists& histos) const {
0712   initOniaToMuMuComponentHistograms(ibook, run, iSetup, histos);
0713   initComponentHists(ibook, run, iSetup, histos, "lead_pi");
0714   initComponentHists(ibook, run, iSetup, histos, "soft_pi");
0715 }
0716 
0717 reco::Vertex const* HeavyFlavorDQMAnalyzer::fillDecayHistograms(DecayHists const& histos,
0718                                                                 pat::CompositeCandidate const& cand,
0719                                                                 reco::VertexCollection const& pvs) const {
0720   //  if (not cand.hasUserData("fitMomentum")) {
0721   //    return -2;
0722   //  }
0723   //  auto mass = cand.userFloat("fitMass");
0724   //  auto& momentum = *cand.userData<GlobalVector>("fitMomentum");
0725   if (not allTracksAvailable(cand)) {
0726     return nullptr;
0727   }
0728 
0729   auto svtx = cand.userData<reco::Vertex>("vertex");
0730   if (not svtx->isValid()) {
0731     return nullptr;
0732   }
0733 
0734   float mass = cand.mass();
0735   reco::Candidate::Vector momentum = cand.momentum();
0736   if (cand.hasUserData("fitMomentum")) {
0737     mass = cand.userFloat("fitMass");
0738     momentum = *cand.userData<GlobalVector>("fitMomentum");
0739   }
0740 
0741   auto pvtx = std::min_element(pvs.begin(), pvs.end(), [svtx](reco::Vertex const& pv1, reco::Vertex const& pv2) {
0742     return abs(pv1.z() - svtx->z()) < abs(pv2.z() - svtx->z());
0743   });
0744   if (pvtx == pvs.end()) {
0745     return nullptr;
0746   }
0747 
0748   VertexDistanceXY vdistXY;
0749   Measurement1D distXY = vdistXY.distance(*svtx, *pvtx);
0750 
0751   auto pvtPos = pvtx->position();
0752   auto svtPos = svtx->position();
0753 
0754   math::XYZVector displVect2D(svtPos.x() - pvtPos.x(), svtPos.y() - pvtPos.y(), 0);
0755   auto cosAlpha = displVect2D.Dot(momentum) / (displVect2D.Rho() * momentum.rho());
0756 
0757   auto ct = distXY.value() * cosAlpha * mass / momentum.rho();
0758 
0759   histos.h_pointing->Fill(cosAlpha);
0760 
0761   histos.h_mass->Fill(mass);
0762 
0763   histos.h_pt->Fill(momentum.rho());
0764   histos.h_eta->Fill(momentum.eta());
0765   histos.h_phi->Fill(momentum.phi());
0766 
0767   histos.h_ct->Fill(ct);
0768 
0769   histos.h_displ2D->Fill(distXY.value());
0770   histos.h_sign2D->Fill(distXY.significance());
0771 
0772   // FIXME workaround for tracks with non pos-def cov. matrix
0773   if (svtx->chi2() >= 0) {
0774     histos.h_vertNormChi2->Fill(svtx->chi2() / svtx->ndof());
0775     histos.h_vertProb->Fill(ChiSquaredProbability(svtx->chi2(), svtx->ndof()));
0776   }
0777 
0778   return &*pvtx;
0779 }
0780 
0781 int HeavyFlavorDQMAnalyzer::fillOniaToMuMuComponents(DecayHists const& histos,
0782                                                      pat::CompositeCandidate const& cand,
0783                                                      reco::BeamSpot const* bs,
0784                                                      reco::Vertex const* pv,
0785                                                      int startPosition) const {
0786   startPosition = fillComponentHistogramsLeadSoft(histos, cand, "MuPos", "MuNeg", bs, pv, startPosition);
0787 
0788   return startPosition;
0789 }
0790 
0791 int HeavyFlavorDQMAnalyzer::fillKx0ToKPiComponents(DecayHists const& histos,
0792                                                    pat::CompositeCandidate const& cand,
0793                                                    reco::BeamSpot const* bs,
0794                                                    reco::Vertex const* pv,
0795                                                    int startPosition) const {
0796   startPosition = fillComponentHistogramsSinglePart(histos, cand, "Kaon", bs, pv, startPosition);
0797   startPosition = fillComponentHistogramsSinglePart(histos, cand, "Pion", bs, pv, startPosition);
0798 
0799   return startPosition;
0800 }
0801 
0802 int HeavyFlavorDQMAnalyzer::fillPhiToKKComponents(DecayHists const& histos,
0803                                                   pat::CompositeCandidate const& cand,
0804                                                   reco::BeamSpot const* bs,
0805                                                   reco::Vertex const* pv,
0806                                                   int startPosition) const {
0807   startPosition = fillComponentHistogramsLeadSoft(histos, cand, "KPos", "KNeg", bs, pv, startPosition);
0808 
0809   return startPosition;
0810 }
0811 
0812 int HeavyFlavorDQMAnalyzer::fillK0sToPiPiComponents(DecayHists const& histos,
0813                                                     pat::CompositeCandidate const& cand,
0814                                                     reco::BeamSpot const* bs,
0815                                                     reco::Vertex const* pv,
0816                                                     int startPosition) const {
0817   startPosition = fillComponentHistogramsLeadSoft(histos, cand, "PionPos", "PionNeg", bs, pv, startPosition);
0818 
0819   return startPosition;
0820 }
0821 
0822 int HeavyFlavorDQMAnalyzer::fillLambda0ToPPiComponents(DecayHists const& histos,
0823                                                        pat::CompositeCandidate const& cand,
0824                                                        reco::BeamSpot const* bs,
0825                                                        reco::Vertex const* pv,
0826                                                        int startPosition) const {
0827   startPosition = fillComponentHistogramsSinglePart(histos, cand, "Proton", bs, pv, startPosition);
0828   startPosition = fillComponentHistogramsSinglePart(histos, cand, "Pion", bs, pv, startPosition);
0829 
0830   return startPosition;
0831 }
0832 
0833 int HeavyFlavorDQMAnalyzer::fillBuToJPsiKComponents(DecayHists const& histos,
0834                                                     pat::CompositeCandidate const& cand,
0835                                                     reco::BeamSpot const* bs,
0836                                                     reco::Vertex const* pv,
0837                                                     int startPosition) const {
0838   startPosition = fillOniaToMuMuComponents(
0839       histos, **cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi"), bs, pv, startPosition);
0840   startPosition = fillComponentHistogramsSinglePart(histos, cand, "Kaon", bs, pv, startPosition);
0841 
0842   return startPosition;
0843 }
0844 
0845 int HeavyFlavorDQMAnalyzer::fillBuToPsi2SKComponents(DecayHists const& histos,
0846                                                      pat::CompositeCandidate const& cand,
0847                                                      reco::BeamSpot const* bs,
0848                                                      reco::Vertex const* pv,
0849                                                      int startPosition) const {
0850   startPosition = fillPsi2SToJPsiPiPiComponents(
0851       histos, **cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToPsi2S"), bs, pv, startPosition);
0852   startPosition = fillComponentHistogramsSinglePart(histos, cand, "Kaon", bs, pv, startPosition);
0853 
0854   return startPosition;
0855 }
0856 
0857 int HeavyFlavorDQMAnalyzer::fillBdToJPsiKx0Components(DecayHists const& histos,
0858                                                       pat::CompositeCandidate const& cand,
0859                                                       reco::BeamSpot const* bs,
0860                                                       reco::Vertex const* pv,
0861                                                       int startPosition) const {
0862   startPosition = fillOniaToMuMuComponents(
0863       histos, **cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi"), bs, pv, startPosition);
0864   startPosition = fillKx0ToKPiComponents(
0865       histos, **cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToKx0"), bs, pv, startPosition);
0866 
0867   return startPosition;
0868 }
0869 
0870 int HeavyFlavorDQMAnalyzer::fillBsToJPsiPhiComponents(DecayHists const& histos,
0871                                                       pat::CompositeCandidate const& cand,
0872                                                       reco::BeamSpot const* bs,
0873                                                       reco::Vertex const* pv,
0874                                                       int startPosition) const {
0875   startPosition = fillOniaToMuMuComponents(
0876       histos, **cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi"), bs, pv, startPosition);
0877   startPosition = fillPhiToKKComponents(
0878       histos, **cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToPhi"), bs, pv, startPosition);
0879 
0880   return startPosition;
0881 }
0882 
0883 int HeavyFlavorDQMAnalyzer::fillBdToJPsiK0sComponents(DecayHists const& histos,
0884                                                       pat::CompositeCandidate const& cand,
0885                                                       reco::BeamSpot const* bs,
0886                                                       reco::Vertex const* pv,
0887                                                       int startPosition) const {
0888   startPosition = fillOniaToMuMuComponents(
0889       histos, **cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi"), bs, pv, startPosition);
0890   startPosition = fillK0sToPiPiComponents(
0891       histos, **cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToK0s"), bs, pv, startPosition);
0892 
0893   return startPosition;
0894 }
0895 
0896 int HeavyFlavorDQMAnalyzer::fillBcToJPsiPiComponents(DecayHists const& histos,
0897                                                      pat::CompositeCandidate const& cand,
0898                                                      reco::BeamSpot const* bs,
0899                                                      reco::Vertex const* pv,
0900                                                      int startPosition) const {
0901   startPosition = fillOniaToMuMuComponents(
0902       histos, **cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi"), bs, pv, startPosition);
0903   startPosition = fillComponentHistogramsSinglePart(histos, cand, "Pion", bs, pv, startPosition);
0904 
0905   return startPosition;
0906 }
0907 
0908 int HeavyFlavorDQMAnalyzer::fillLambdaBToJPsiLambda0Components(DecayHists const& histos,
0909                                                                pat::CompositeCandidate const& cand,
0910                                                                reco::BeamSpot const* bs,
0911                                                                reco::Vertex const* pv,
0912                                                                int startPosition) const {
0913   startPosition = fillOniaToMuMuComponents(
0914       histos, **cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi"), bs, pv, startPosition);
0915   startPosition = fillLambda0ToPPiComponents(
0916       histos, **cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToLambda0"), bs, pv, startPosition);
0917 
0918   return startPosition;
0919 }
0920 
0921 int HeavyFlavorDQMAnalyzer::fillPsi2SToJPsiPiPiComponents(DecayHists const& histos,
0922                                                           pat::CompositeCandidate const& cand,
0923                                                           reco::BeamSpot const* bs,
0924                                                           reco::Vertex const* pv,
0925                                                           int startPosition) const {
0926   startPosition = fillOniaToMuMuComponents(
0927       histos, **cand.userData<edm::Ref<pat::CompositeCandidateCollection>>("refToJPsi"), bs, pv, startPosition);
0928   startPosition = fillComponentHistogramsLeadSoft(histos, cand, "PionPos", "PionNeg", bs, pv, startPosition);
0929   return startPosition;
0930 }
0931 
0932 void HeavyFlavorDQMAnalyzer::fillComponentHistograms(ComponentHists const& histos,
0933                                                      reco::Track const& component,
0934                                                      reco::BeamSpot const* bs,
0935                                                      reco::Vertex const* pv) const {
0936   histos.h_pt->Fill(component.pt());
0937   histos.h_eta->Fill(component.eta());
0938   histos.h_phi->Fill(component.phi());
0939 
0940   math::XYZPoint zero(0, 0, 0);
0941   math::Error<3>::type zeroCov;  // needed for dxyError
0942   if (bs) {
0943     histos.h_dxy->Fill(component.dxy(*bs));
0944     histos.h_exy->Fill(component.dxyError(*bs));
0945   } else {
0946     histos.h_dxy->Fill(component.dxy(zero));
0947     histos.h_exy->Fill(component.dxyError(zero, zeroCov));
0948   }
0949   if (pv) {
0950     histos.h_dz->Fill(component.dz(pv->position()));
0951   } else {
0952     histos.h_dz->Fill(component.dz(zero));
0953   }
0954   histos.h_ez->Fill(component.dzError());
0955 
0956   histos.h_chi2->Fill(component.chi2() / component.ndof());
0957 }
0958 
0959 bool HeavyFlavorDQMAnalyzer::allTracksAvailable(pat::CompositeCandidate const& cand) const {
0960   for (auto&& name : cand.roles()) {
0961     auto track = getDaughterTrack(cand, name, false);
0962     if (not track) {
0963       return false;
0964     }
0965   }
0966   return true;
0967 }
0968 
0969 const reco::Track* HeavyFlavorDQMAnalyzer::getDaughterTrack(pat::CompositeCandidate const& cand,
0970                                                             std::string const& name,
0971                                                             bool throwOnMissing) const {
0972   auto daugh = cand.daughter(name);
0973   auto trackModeLabel = "trackMode_" + name;
0974   auto trackMode = cand.userData<std::string>(trackModeLabel);
0975   if (!trackMode or trackMode->empty()) {
0976     if (throwOnMissing) {
0977       throw cms::Exception("TrackNotFound") << "Could not determine track mode from candidate with name " << name
0978                                             << " with label " << trackModeLabel << std::endl;
0979     }
0980     return nullptr;
0981   }
0982 
0983   auto track = BPHTrackReference::getTrack(*daugh, trackMode->c_str());
0984 
0985   if (throwOnMissing and not track) {
0986     throw cms::Exception("TrackNotFound") << "BPHTrackReference could not extract a track as type " << trackMode
0987                                           << " from candidate with name " << name << std::endl;
0988   }
0989 
0990   return track;
0991 }
0992 
0993 int HeavyFlavorDQMAnalyzer::fillComponentHistogramsSinglePart(DecayHists const& histos,
0994                                                               pat::CompositeCandidate const& cand,
0995                                                               std::string const& name,
0996                                                               reco::BeamSpot const* bs,
0997                                                               reco::Vertex const* pv,
0998                                                               int startPosition) const {
0999   fillComponentHistograms(histos.decayComponents[startPosition], *getDaughterTrack(cand, name), bs, pv);
1000 
1001   return startPosition + 1;
1002 }
1003 
1004 int HeavyFlavorDQMAnalyzer::fillComponentHistogramsLeadSoft(DecayHists const& histos,
1005                                                             pat::CompositeCandidate const& cand,
1006                                                             std::string const& name1,
1007                                                             std::string const& name2,
1008                                                             reco::BeamSpot const* bs,
1009                                                             reco::Vertex const* pv,
1010                                                             int startPosition) const {
1011   auto daughSoft = getDaughterTrack(cand, name1);
1012   auto daughLead = getDaughterTrack(cand, name2);
1013 
1014   if (daughLead->pt() < daughSoft->pt()) {
1015     std::swap(daughLead, daughSoft);
1016   }
1017 
1018   fillComponentHistograms(histos.decayComponents[startPosition], *daughLead, bs, pv);
1019   fillComponentHistograms(histos.decayComponents[startPosition + 1], *daughSoft, bs, pv);
1020 
1021   return startPosition + 2;
1022 }
1023 
1024 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
1025 void HeavyFlavorDQMAnalyzer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
1026   edm::ParameterSetDescription desc;
1027 
1028   desc.add<std::string>("folder", "Physics/HeavyFlavor");
1029 
1030   desc.add<edm::InputTag>("pvCollection");
1031   desc.add<edm::InputTag>("beamSpot");
1032 
1033   desc.addOptional<edm::InputTag>("OniaToMuMuCands");
1034   desc.addOptional<edm::InputTag>("Kx0ToKPiCands");
1035   desc.addOptional<edm::InputTag>("PhiToKKCands");
1036   desc.addOptional<edm::InputTag>("BuToJPsiKCands");
1037   desc.addOptional<edm::InputTag>("BuToPsi2SKCands");
1038   desc.addOptional<edm::InputTag>("BdToJPsiKx0Cands");
1039   desc.addOptional<edm::InputTag>("BsToJPsiPhiCands");
1040   desc.addOptional<edm::InputTag>("K0sToPiPiCands");
1041   desc.addOptional<edm::InputTag>("Lambda0ToPPiCands");
1042   desc.addOptional<edm::InputTag>("BdToJPsiK0sCands");
1043   desc.addOptional<edm::InputTag>("LambdaBToJPsiLambda0Cands");
1044   desc.addOptional<edm::InputTag>("BcToJPsiPiCands");
1045   desc.addOptional<edm::InputTag>("Psi2SToJPsiPiPiCands");
1046 
1047   descriptions.add("HeavyFlavorDQMAnalyzer", desc);
1048 }
1049 
1050 // define this as a plug-in
1051 DEFINE_FWK_MODULE(HeavyFlavorDQMAnalyzer);