ComponentHists

DecayHists

HeavyFlavorDQMAnalyzer

Histograms_HeavyFlavorDQMAnalyzer

Macros

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283
#ifndef DQMOffline_Physics_HeavyFlavorDQMAnalyzer_h
#define DQMOffline_Physics_HeavyFlavorDQMAnalyzer_h

// -*- C++ -*-
//
// Package:    DQMOffline/HeavyFlavorDQMAnalyzer
// Class:      HeavyFlavorDQMAnalyzer
//
/**\class HeavyFlavorDQMAnalyzer HeavyFlavorDQMAnalyzer.cc DQMOffline/HeavyFlavorDQMAnalyzer/plugins/HeavyFlavorDQMAnalyzer.cc

 Description: [one line class summary]

 Implementation:
     [Notes on implementation]
*/
//
// Original Author:  Enrico Lusiani
//         Created:  Mon, 22 Nov 2021 14:36:39 GMT
//
//

#include <string>

// user include files
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "DQMServices/Core/interface/DQMGlobalEDAnalyzer.h"

#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"

#include "FWCore/ParameterSet/interface/ParameterSet.h"

#include "DataFormats/Common/interface/Ref.h"
#include "DataFormats/PatCandidates/interface/CompositeCandidate.h"
#include "DataFormats/VertexReco/interface/Vertex.h"

#include "CommonTools/Statistics/interface/ChiSquaredProbability.h"
#include "RecoVertex/VertexTools/interface/VertexDistanceXY.h"

#include "HeavyFlavorAnalysis/RecoDecay/interface/BPHTrackReference.h"

//
// class declaration
//

struct ComponentHists {
  dqm::reco::MonitorElement* h_pt;
  dqm::reco::MonitorElement* h_eta;
  dqm::reco::MonitorElement* h_phi;

  dqm::reco::MonitorElement* h_dxy;
  dqm::reco::MonitorElement* h_exy;
  dqm::reco::MonitorElement* h_dz;
  dqm::reco::MonitorElement* h_ez;

  dqm::reco::MonitorElement* h_chi2;
};

struct DecayHists {
  // kinematics
  dqm::reco::MonitorElement* h_mass;
  dqm::reco::MonitorElement* h_pt;
  dqm::reco::MonitorElement* h_eta;
  dqm::reco::MonitorElement* h_phi;

  // position
  dqm::reco::MonitorElement* h_displ2D;
  dqm::reco::MonitorElement* h_displ3D;
  dqm::reco::MonitorElement* h_sign2D;
  dqm::reco::MonitorElement* h_sign3D;

  // ct and pointing angle
  dqm::reco::MonitorElement* h_ct;
  dqm::reco::MonitorElement* h_pointing;

  // quality
  dqm::reco::MonitorElement* h_vertNormChi2;
  dqm::reco::MonitorElement* h_vertProb;

  std::vector<ComponentHists> decayComponents;
};

struct Histograms_HeavyFlavorDQMAnalyzer {
  DecayHists oniaToMuMuPrompt;
  DecayHists oniaToMuMuDispl;
  DecayHists oniaToMuMu;
  DecayHists kx0ToKPiPrompt;
  DecayHists kx0ToKPiDispl;
  DecayHists kx0ToKPi;
  DecayHists phiToKKPrompt;
  DecayHists phiToKKDispl;
  DecayHists phiToKK;
  DecayHists psi2SToJPsiPiPiPrompt;
  DecayHists psi2SToJPsiPiPiDispl;
  DecayHists psi2SToJPsiPiPi;
  DecayHists k0sToPiPi;
  DecayHists lambda0ToPPi;
  DecayHists buToJPsiK;
  DecayHists buToPsi2SK;
  DecayHists bdToJPsiKx0;
  DecayHists bsToJPsiPhi;
  DecayHists bdToJPsiK0s;
  DecayHists bcToJPsiPi;
  DecayHists lambdaBToJPsiLambda0;
};

class HeavyFlavorDQMAnalyzer : public DQMGlobalEDAnalyzer<Histograms_HeavyFlavorDQMAnalyzer> {
public:
  using Histograms = Histograms_HeavyFlavorDQMAnalyzer;

  explicit HeavyFlavorDQMAnalyzer(const edm::ParameterSet&);
  ~HeavyFlavorDQMAnalyzer() override;

  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

private:
  void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&, Histograms&) const override;

  void dqmAnalyze(edm::Event const&, edm::EventSetup const&, Histograms const&) const override;

  void bookDecayHists(DQMStore::IBooker&,
                      edm::Run const&,
                      edm::EventSetup const&,
                      DecayHists&,
                      std::string const&,
                      std::string const&,
                      int,
                      float,
                      float,
                      float distanceScaleFactor = 1.) const;
  void initComponentHists(DQMStore::IBooker&,
                          edm::Run const&,
                          edm::EventSetup const&,
                          DecayHists&,
                          TString const&) const;  // TString for the IBooker interface

  void initOniaToMuMuComponentHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&, DecayHists&) const;
  void initKx0ToKPiComponentHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&, DecayHists&) const;
  void initPhiToKKComponentHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&, DecayHists&) const;
  void initK0sToPiPiComponentHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&, DecayHists&) const;
  void initLambda0ToPPiComponentHistograms(DQMStore::IBooker&,
                                           edm::Run const&,
                                           edm::EventSetup const&,
                                           DecayHists&) const;
  void initBuToJPsiKComponentHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&, DecayHists&) const;
  void initBuToPsi2SKComponentHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&, DecayHists&) const;
  void initBdToJPsiKx0ComponentHistograms(DQMStore::IBooker&,
                                          edm::Run const&,
                                          edm::EventSetup const&,
                                          DecayHists&) const;
  void initBsToJPsiPhiComponentHistograms(DQMStore::IBooker&,
                                          edm::Run const&,
                                          edm::EventSetup const&,
                                          DecayHists&) const;
  void initBdToJPsiK0sComponentHistograms(DQMStore::IBooker&,
                                          edm::Run const&,
                                          edm::EventSetup const&,
                                          DecayHists&) const;
  void initBcToJPsiPiComponentHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&, DecayHists&) const;
  void initLambdaBToJPsiLambda0ComponentHistograms(DQMStore::IBooker&,
                                                   edm::Run const&,
                                                   edm::EventSetup const&,
                                                   DecayHists&) const;
  void initPsi2SToJPsiPiPiComponentHistograms(DQMStore::IBooker&,
                                              edm::Run const&,
                                              edm::EventSetup const&,
                                              DecayHists&) const;

  reco::Vertex const* fillDecayHistograms(DecayHists const&,
                                          pat::CompositeCandidate const& cand,
                                          reco::VertexCollection const& pvs) const;
  void fillComponentHistograms(ComponentHists const& histos,
                               reco::Track const& component,
                               reco::BeamSpot const* bs,
                               reco::Vertex const* pv) const;

  int fillComponentHistogramsSinglePart(DecayHists const&,
                                        pat::CompositeCandidate const& cand,
                                        std::string const& name,
                                        reco::BeamSpot const* bs,
                                        reco::Vertex const* pv,
                                        int startPosition = 0) const;
  int fillComponentHistogramsLeadSoft(DecayHists const&,
                                      pat::CompositeCandidate const& cand,
                                      std::string const& name1,
                                      std::string const& name2,
                                      reco::BeamSpot const* bs,
                                      reco::Vertex const* pv,
                                      int startPosition = 0) const;

  const reco::Track* getDaughterTrack(pat::CompositeCandidate const& cand,
                                      std::string const& name,
                                      bool throwOnMissing = true) const;
  bool allTracksAvailable(pat::CompositeCandidate const& cand) const;

  int fillOniaToMuMuComponents(DecayHists const& histos,
                               pat::CompositeCandidate const& cand,
                               reco::BeamSpot const* bs,
                               reco::Vertex const* pv,
                               int startPosition = 0) const;
  int fillKx0ToKPiComponents(DecayHists const& histos,
                             pat::CompositeCandidate const& cand,
                             reco::BeamSpot const* bs,
                             reco::Vertex const* pv,
                             int startPosition = 0) const;
  int fillPhiToKKComponents(DecayHists const& histos,
                            pat::CompositeCandidate const& cand,
                            reco::BeamSpot const* bs,
                            reco::Vertex const* pv,
                            int startPosition = 0) const;
  int fillK0sToPiPiComponents(DecayHists const& histos,
                              pat::CompositeCandidate const& cand,
                              reco::BeamSpot const* bs,
                              reco::Vertex const* pv,
                              int startPosition = 0) const;
  int fillLambda0ToPPiComponents(DecayHists const& histos,
                                 pat::CompositeCandidate const& cand,
                                 reco::BeamSpot const* bs,
                                 reco::Vertex const* pv,
                                 int startPosition = 0) const;
  int fillBuToJPsiKComponents(DecayHists const& histos,
                              pat::CompositeCandidate const& cand,
                              reco::BeamSpot const* bs,
                              reco::Vertex const* pv,
                              int startPosition = 0) const;
  int fillBuToPsi2SKComponents(DecayHists const& histos,
                               pat::CompositeCandidate const& cand,
                               reco::BeamSpot const* bs,
                               reco::Vertex const* pv,
                               int startPosition = 0) const;
  int fillBdToJPsiKx0Components(DecayHists const& histos,
                                pat::CompositeCandidate const& cand,
                                reco::BeamSpot const* bs,
                                reco::Vertex const* pv,
                                int startPosition = 0) const;
  int fillBsToJPsiPhiComponents(DecayHists const& histos,
                                pat::CompositeCandidate const& cand,
                                reco::BeamSpot const* bs,
                                reco::Vertex const* pv,
                                int startPosition = 0) const;
  int fillBdToJPsiK0sComponents(DecayHists const& histos,
                                pat::CompositeCandidate const& cand,
                                reco::BeamSpot const* bs,
                                reco::Vertex const* pv,
                                int startPosition = 0) const;
  int fillBcToJPsiPiComponents(DecayHists const& histos,
                               pat::CompositeCandidate const& cand,
                               reco::BeamSpot const* bs,
                               reco::Vertex const* pv,
                               int startPosition = 0) const;
  int fillLambdaBToJPsiLambda0Components(DecayHists const& histos,
                                         pat::CompositeCandidate const& cand,
                                         reco::BeamSpot const* bs,
                                         reco::Vertex const* pv,
                                         int startPosition = 0) const;
  int fillPsi2SToJPsiPiPiComponents(DecayHists const& histos,
                                    pat::CompositeCandidate const& cand,
                                    reco::BeamSpot const* bs,
                                    reco::Vertex const* pv,
                                    int startPosition = 0) const;

  // ------------ member data ------------
  std::string folder_;

  edm::EDGetTokenT<reco::VertexCollection> pvCollectionToken;
  edm::EDGetTokenT<reco::BeamSpot> beamSpotToken;

  edm::EDGetTokenT<pat::CompositeCandidateCollection> oniaToMuMuCandsToken;
  edm::EDGetTokenT<pat::CompositeCandidateCollection> kx0ToKPiCandsToken;
  edm::EDGetTokenT<pat::CompositeCandidateCollection> phiToKKCandsToken;
  edm::EDGetTokenT<pat::CompositeCandidateCollection> buToJPsiKCandsToken;
  edm::EDGetTokenT<pat::CompositeCandidateCollection> buToPsi2SKCandsToken;
  edm::EDGetTokenT<pat::CompositeCandidateCollection> bdToJPsiKx0CandsToken;
  edm::EDGetTokenT<pat::CompositeCandidateCollection> bsToJPsiPhiCandsToken;
  edm::EDGetTokenT<pat::CompositeCandidateCollection> k0sToPiPiCandsToken;
  edm::EDGetTokenT<pat::CompositeCandidateCollection> lambda0ToPPiCandsToken;
  edm::EDGetTokenT<pat::CompositeCandidateCollection> bdToJPsiK0sCandsToken;
  edm::EDGetTokenT<pat::CompositeCandidateCollection> lambdaBToJPsiLambda0CandsToken;
  edm::EDGetTokenT<pat::CompositeCandidateCollection> bcToJPsiPiCandsToken;
  edm::EDGetTokenT<pat::CompositeCandidateCollection> psi2SToJPsiPiPiCandsToken;
};

#endif