CMS 3D CMS Logo

List of all members | Public Member Functions | Static Public Member Functions | Private Attributes
DiMuonMassBiasMonitor Class Reference

#include <DiMuonMassBiasMonitor.h>

Inheritance diagram for DiMuonMassBiasMonitor:
DQMEDAnalyzer edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &) override
 
void bookComponentHists (DQMStore::IBooker &, DecayHists &, TString const &, float distanceScaleFactor=1.) const
 
void bookDecayComponentHistograms (DQMStore::IBooker &ibook, DecayHists &histos) const
 
void bookDecayHists (DQMStore::IBooker &, DecayHists &, std::string const &, std::string const &, int, float, float, float distanceScaleFactor=1.) const
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
 DiMuonMassBiasMonitor (const edm::ParameterSet &)
 
void fillComponentHistograms (ComponentHists const &histos, const reco::Track *const &component, reco::BeamSpot const *bs, reco::Vertex const *pv) const
 
reco::Vertex const * fillDecayHistograms (DecayHists const &, std::vector< const reco::Track *> const &tracks, const reco::VertexCollection *const &pvs, const edm::EventSetup &) const
 
 ~DiMuonMassBiasMonitor () override=default
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &event, edm::EventSetup const &setup) final
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
void beginStream (edm::StreamID id) final
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 DQMEDAnalyzer ()
 
void endLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual bool getCanSaveByLumi ()
 
- Public Member Functions inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
 EDProducer ()=default
 
 EDProducer (const EDProducer &)=delete
 
bool hasAbilityToProduceInBeginLumis () const final
 
bool hasAbilityToProduceInBeginProcessBlocks () const final
 
bool hasAbilityToProduceInBeginRuns () const final
 
bool hasAbilityToProduceInEndLumis () const final
 
bool hasAbilityToProduceInEndProcessBlocks () const final
 
bool hasAbilityToProduceInEndRuns () const final
 
const EDProduceroperator= (const EDProducer &)=delete
 

Static Public Member Functions

static void fillDescriptions (edm::ConfigurationDescriptions &descriptions)
 
- Static Public Member Functions inherited from DQMEDAnalyzer
static void globalEndJob (DQMEDAnalyzerGlobalCache const *)
 
static void globalEndLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup, LuminosityBlockContext const *context)
 
static void globalEndRunProduce (edm::Run &run, edm::EventSetup const &setup, RunContext const *context)
 
static std::unique_ptr< DQMEDAnalyzerGlobalCacheinitializeGlobalCache (edm::ParameterSet const &)
 

Private Attributes

const edm::EDGetTokenT< reco::BeamSpotbeamSpotToken_
 
const std::string decayMotherName_
 
edm::ParameterSet DiMuMassConfiguration_
 
const double distanceScaleFactor_
 
DecayHists histosZmm
 
const std::string MEFolderName_
 
const edm::EDGetTokenT< reco::TrackCollectiontracksToken_
 
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecordttbESToken_
 
const edm::EDGetTokenT< reco::VertexCollectionvertexToken_
 
DiLepPlotHelp::PlotsVsKinematics ZMassPlots = DiLepPlotHelp::PlotsVsKinematics(DiLepPlotHelp::MM)
 

Additional Inherited Members

- Public Types inherited from DQMEDAnalyzer
typedef dqm::reco::DQMStore DQMStore
 
typedef dqm::reco::MonitorElement MonitorElement
 
- Public Types inherited from edm::stream::EDProducer< edm::GlobalCache< DQMEDAnalyzerGlobalCache >, edm::EndRunProducer, edm::EndLuminosityBlockProducer, edm::Accumulator >
using CacheTypes = CacheContexts< T... >
 
using GlobalCache = typename CacheTypes::GlobalCache
 
using HasAbility = AbilityChecker< T... >
 
using InputProcessBlockCache = typename CacheTypes::InputProcessBlockCache
 
using LuminosityBlockCache = typename CacheTypes::LuminosityBlockCache
 
using LuminosityBlockContext = LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCache >
 
using LuminosityBlockSummaryCache = typename CacheTypes::LuminosityBlockSummaryCache
 
using RunCache = typename CacheTypes::RunCache
 
using RunContext = RunContextT< RunCache, GlobalCache >
 
using RunSummaryCache = typename CacheTypes::RunSummaryCache
 
- Protected Member Functions inherited from DQMEDAnalyzer
uint64_t meId () const
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 
unsigned int streamId_
 

Detailed Description

DQM/TrackerMonitorTrack/src/DiMuonMassBiasMonitor.cc Monitoring DiMuon mass bias

Definition at line 69 of file DiMuonMassBiasMonitor.h.

Constructor & Destructor Documentation

◆ DiMuonMassBiasMonitor()

DiMuonMassBiasMonitor::DiMuonMassBiasMonitor ( const edm::ParameterSet iConfig)
explicit

Definition at line 29 of file DiMuonMassBiasMonitor.cc.

30  : ttbESToken_(esConsumes(edm::ESInputTag("", "TransientTrackBuilder"))),
31  tracksToken_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("muonTracks"))),
32  vertexToken_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("vertices"))),
33  beamSpotToken_(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("beamSpot"))),
34  MEFolderName_(iConfig.getParameter<std::string>("FolderName")),
35  decayMotherName_(iConfig.getParameter<std::string>("decayMotherName")),
36  distanceScaleFactor_(iConfig.getParameter<double>("distanceScaleFactor")),
37  DiMuMassConfiguration_(iConfig.getParameter<edm::ParameterSet>("DiMuMassConfig")) {}
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
edm::ParameterSet DiMuMassConfiguration_
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > ttbESToken_
const edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
const edm::EDGetTokenT< reco::VertexCollection > vertexToken_
const edm::EDGetTokenT< reco::TrackCollection > tracksToken_
const std::string decayMotherName_
const std::string MEFolderName_

◆ ~DiMuonMassBiasMonitor()

DiMuonMassBiasMonitor::~DiMuonMassBiasMonitor ( )
overridedefault

Member Function Documentation

◆ analyze()

void DiMuonMassBiasMonitor::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overridevirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 57 of file DiMuonMassBiasMonitor.cc.

References beamSpotToken_, cms::cuda::bs, DecayHists::decayComponents, fillComponentHistograms(), fillDecayHistograms(), DiLepPlotHelp::PlotsVsKinematics::fillPlots(), histosZmm, iEvent, mumass2, mathSSE::sqrt(), FrontierCondition_GT_autoExpress_cfi::t0, RandomServiceHelper::t1, tracksToken_, vertexToken_, AlignmentTracksFromVertexSelector_cfi::vertices, and ZMassPlots.

57  {
58  std::vector<const reco::Track*> myTracks;
59  const auto trackHandle = iEvent.getHandle(tracksToken_);
60  if (!trackHandle.isValid()) {
61  edm::LogError("DiMuonMassBiasMonitor") << "invalid track collection encountered!";
62  return;
63  }
64 
65  for (const auto& muonTrk : *trackHandle) {
66  myTracks.emplace_back(&muonTrk);
67  }
68 
69  if (myTracks.size() != 2) {
70  edm::LogWarning("DiMuonMassBiasMonitor") << "There are not enough tracks to monitor!";
71  return;
72  }
73 
74  const auto& t1 = myTracks[1]->momentum();
75  const auto& t0 = myTracks[0]->momentum();
76  const auto& ditrack = t1 + t0;
77 
78  const auto& tplus = myTracks[0]->charge() > 0 ? myTracks[0] : myTracks[1];
79  const auto& tminus = myTracks[0]->charge() < 0 ? myTracks[0] : myTracks[1];
80 
81  TLorentzVector p4_tplus(tplus->px(), tplus->py(), tplus->pz(), sqrt((tplus->p() * tplus->p()) + mumass2));
82  TLorentzVector p4_tminus(tminus->px(), tminus->py(), tminus->pz(), sqrt((tminus->p() * tminus->p()) + mumass2));
83 
84  const auto& Zp4 = p4_tplus + p4_tminus;
85  float track_invMass = Zp4.M();
86 
87  // creat the pair of TLorentVectors used to make the plos
88  std::pair<TLorentzVector, TLorentzVector> tktk_p4 = std::make_pair(p4_tplus, p4_tminus);
89 
90  // fill the z->mm mass plots
91  ZMassPlots.fillPlots(track_invMass, tktk_p4);
92 
93  math::XYZPoint ZpT(ditrack.x(), ditrack.y(), 0);
94  math::XYZPoint Zp(ditrack.x(), ditrack.y(), ditrack.z());
95 
96  // get collection of reconstructed vertices from event
97  const auto& vertexHandle = iEvent.getHandle(vertexToken_);
98  if (!vertexHandle.isValid()) {
99  edm::LogError("DiMuonMassBiasMonitor") << "invalid vertex collection encountered!";
100  return;
101  }
102 
103  // get the vertices from the event
104  const auto& vertices = vertexHandle.product();
105 
106  // fill the decay vertex plots
107  auto decayVertex = fillDecayHistograms(histosZmm, myTracks, vertices, iSetup);
108 
109  // get the beamspot from the event
110  const reco::BeamSpot* bs;
111  const auto& beamSpotHandle = iEvent.getHandle(beamSpotToken_);
112  if (!beamSpotHandle.isValid()) {
113  bs = nullptr;
114  } else {
115  bs = beamSpotHandle.product();
116  }
117 
118  // fill the components plots
119  fillComponentHistograms(histosZmm.decayComponents[0], tplus, bs, decayVertex);
120  fillComponentHistograms(histosZmm.decayComponents[1], tminus, bs, decayVertex);
121 }
static constexpr float mumass2
Log< level::Error, false > LogError
reco::Vertex const * fillDecayHistograms(DecayHists const &, std::vector< const reco::Track *> const &tracks, const reco::VertexCollection *const &pvs, const edm::EventSetup &) const
int iEvent
Definition: GenABIO.cc:224
T sqrt(T t)
Definition: SSEVec.h:19
void fillPlots(const float val, const std::pair< TLorentzVector, TLorentzVector > &momenta)
DiLepPlotHelp::PlotsVsKinematics ZMassPlots
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
const edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
const edm::EDGetTokenT< reco::VertexCollection > vertexToken_
Log< level::Warning, false > LogWarning
const edm::EDGetTokenT< reco::TrackCollection > tracksToken_
void fillComponentHistograms(ComponentHists const &histos, const reco::Track *const &component, reco::BeamSpot const *bs, reco::Vertex const *pv) const
std::vector< ComponentHists > decayComponents

◆ bookComponentHists()

void DiMuonMassBiasMonitor::bookComponentHists ( DQMStore::IBooker ibook,
DecayHists histos,
TString const &  componentName,
float  distanceScaleFactor = 1. 
) const

Definition at line 128 of file DiMuonMassBiasMonitor.cc.

References dqm::implementation::IBooker::book1D(), AlCaHLTBitMon_QueryRunRegistry::comp, DiMuonMassBiasMonitor_cfi::distanceScaleFactor, combine::histos, and M_PI.

Referenced by bookDecayComponentHistograms().

131  {
133 
134  comp.h_pt = ibook.book1D(componentName + "_pt", "track momentum ;p_{T} [GeV]", 100, 0., 100.);
135  comp.h_eta = ibook.book1D(componentName + "_eta", "track rapidity;#eta", 100, -2.5, 2.5);
136  comp.h_phi = ibook.book1D(componentName + "_phi", "track azimuth;#phi", 100, -M_PI, M_PI);
137  comp.h_dxy = ibook.book1D(componentName + "_dxyBS",
138  "TIP w.r.t BS;d_{xy}(BS) [cm]",
139  100,
140  -0.5 * distanceScaleFactor,
141  0.5 * distanceScaleFactor);
142  comp.h_exy =
143  ibook.book1D(componentName + "_exy", "Error on TIP ;#sigma(d_{xy}(BS)) [cm]", 100, 0, 0.05 * distanceScaleFactor);
144  comp.h_dz = ibook.book1D(componentName + "_dzPV",
145  "LIP w.r.t PV;d_{z}(PV) [cm]",
146  100,
147  -0.5 * distanceScaleFactor,
148  0.5 * distanceScaleFactor);
149  comp.h_ez =
150  ibook.book1D(componentName + "_ez", "Error on LIP;#sigma(d_{z}(PV)) [cm]", 100, 0, 0.5 * distanceScaleFactor);
151  comp.h_chi2 = ibook.book1D(componentName + "_chi2", ";#chi^{2}", 100, 0, 20);
152 
153  histos.decayComponents.push_back(comp);
154 }
#define M_PI
histos
Definition: combine.py:4
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98

◆ bookDecayComponentHistograms()

void DiMuonMassBiasMonitor::bookDecayComponentHistograms ( DQMStore::IBooker ibook,
DecayHists histos 
) const

Definition at line 123 of file DiMuonMassBiasMonitor.cc.

References bookComponentHists(), and combine::histos.

Referenced by bookHistograms().

123  {
124  bookComponentHists(ibook, histos, "mu_plus", 0.1);
125  bookComponentHists(ibook, histos, "mu_minus", 0.1);
126 }
histos
Definition: combine.py:4
void bookComponentHists(DQMStore::IBooker &, DecayHists &, TString const &, float distanceScaleFactor=1.) const

◆ bookDecayHists()

void DiMuonMassBiasMonitor::bookDecayHists ( DQMStore::IBooker ibook,
DecayHists decayHists,
std::string const &  name,
std::string const &  products,
int  nMassBins,
float  massMin,
float  massMax,
float  distanceScaleFactor = 1. 
) const

Definition at line 156 of file DiMuonMassBiasMonitor.cc.

References dqm::implementation::IBooker::book1D(), DiMuonMassBiasMonitor_cfi::distanceScaleFactor, DecayHists::h_ct, DecayHists::h_displ2D, DecayHists::h_eta, DecayHists::h_mass, DecayHists::h_phi, DecayHists::h_pointing, DecayHists::h_pt, DecayHists::h_sign2D, DecayHists::h_vertNormChi2, DecayHists::h_vertProb, B2GTnPMonitor_cfi::histTitle, M_PI, vertexSelectForHeavyFlavorDQM_cfi::massMax, vertexSelectForHeavyFlavorDQM_cfi::massMin, Skims_PA_cff::name, edm::es::products(), and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by bookHistograms().

163  {
164  std::string histTitle = name + " #rightarrow " + products + ";";
165 
166  decayHists.h_mass = ibook.book1D("h_mass", histTitle + "M(" + products + ") [GeV]", nMassBins, massMin, massMax);
167  decayHists.h_pt = ibook.book1D("h_pt", histTitle + "p_{T} [GeV]", 100, 0.00, 200.0);
168  decayHists.h_eta = ibook.book1D("h_eta", histTitle + "#eta", 100, -5., 5.);
169  decayHists.h_phi = ibook.book1D("h_phi", histTitle + "#varphi [rad]", 100, -M_PI, M_PI);
170  decayHists.h_displ2D =
171  ibook.book1D("h_displ2D", histTitle + "vertex 2D displacement [cm]", 100, 0.00, 0.1 * distanceScaleFactor);
172  decayHists.h_sign2D =
173  ibook.book1D("h_sign2D", histTitle + "vertex 2D displ. significance", 100, 0.00, 100.0 * distanceScaleFactor);
174  decayHists.h_ct = ibook.book1D("h_ct", histTitle + "c#tau [cm]", 100, 0.00, 0.4 * distanceScaleFactor);
175  decayHists.h_pointing = ibook.book1D("h_pointing", histTitle + "cos( 2D pointing angle )", 100, -1, 1);
176  decayHists.h_vertNormChi2 = ibook.book1D("h_vertNormChi2", histTitle + "vertex #chi^{2}/ndof", 100, 0.00, 10);
177  decayHists.h_vertProb = ibook.book1D("h_vertProb", histTitle + "vertex prob.", 100, 0.00, 1.0);
178 }
dqm::reco::MonitorElement * h_displ2D
dqm::reco::MonitorElement * h_phi
ESProducts< std::remove_reference_t< TArgs >... > products(TArgs &&... args)
Definition: ESProducts.h:128
dqm::reco::MonitorElement * h_vertNormChi2
dqm::reco::MonitorElement * h_vertProb
dqm::reco::MonitorElement * h_pointing
dqm::reco::MonitorElement * h_ct
dqm::reco::MonitorElement * h_mass
#define M_PI
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
dqm::reco::MonitorElement * h_eta
dqm::reco::MonitorElement * h_sign2D
dqm::reco::MonitorElement * h_pt

◆ bookHistograms()

void DiMuonMassBiasMonitor::bookHistograms ( DQMStore::IBooker iBooker,
edm::Run const &  ,
edm::EventSetup const &   
)
overridevirtual

Implements DQMEDAnalyzer.

Definition at line 39 of file DiMuonMassBiasMonitor.cc.

References bookDecayComponentHistograms(), bookDecayHists(), DiLepPlotHelp::PlotsVsKinematics::bookFromPSet(), decayMotherName_, DiMuMassConfiguration_, distanceScaleFactor_, edm::ParameterSet::getParameter(), histosZmm, MEFolderName_, dqm::implementation::NavigatorBase::setCurrentFolder(), and ZMassPlots.

39  {
40  iBooker.setCurrentFolder(MEFolderName_ + "/DiMuonMassBiasMonitor/MassBias");
42 
43  iBooker.setCurrentFolder(MEFolderName_ + "/DiMuonMassBiasMonitor");
44 
45  // retrieve the mass bins
46  const auto& mass_bins = DiMuMassConfiguration_.getParameter<int32_t>("NyBins");
47  const auto& mass_min = DiMuMassConfiguration_.getParameter<double>("ymin");
48  const auto& mass_max = DiMuMassConfiguration_.getParameter<double>("ymax");
49 
51  iBooker, histosZmm, decayMotherName_, "#mu^{+}#mu^{-}", mass_bins, mass_min, mass_max, distanceScaleFactor_);
52 
53  iBooker.setCurrentFolder(MEFolderName_ + "/DiMuonMassBiasMonitor/components");
55 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
void bookDecayHists(DQMStore::IBooker &, DecayHists &, std::string const &, std::string const &, int, float, float, float distanceScaleFactor=1.) const
edm::ParameterSet DiMuMassConfiguration_
void bookDecayComponentHistograms(DQMStore::IBooker &ibook, DecayHists &histos) const
DiLepPlotHelp::PlotsVsKinematics ZMassPlots
void bookFromPSet(dqm::reco::DQMStore::IBooker &iBooker, const edm::ParameterSet &hpar)
const std::string decayMotherName_
const std::string MEFolderName_

◆ fillComponentHistograms()

void DiMuonMassBiasMonitor::fillComponentHistograms ( ComponentHists const &  histos,
const reco::Track *const &  component,
reco::BeamSpot const *  bs,
reco::Vertex const *  pv 
) const

Definition at line 256 of file DiMuonMassBiasMonitor.cc.

References cms::cuda::bs, reco::TrackBase::chi2(), reco::TrackBase::dxy(), reco::TrackBase::dxyError(), reco::TrackBase::dz(), reco::TrackBase::dzError(), reco::TrackBase::eta(), combine::histos, reco::TrackBase::ndof(), reco::TrackBase::phi(), reco::TrackBase::pt(), MetAnalyzer::pv(), and SiPixelPI::zero.

Referenced by analyze().

259  {
260  histos.h_pt->Fill(component->pt());
261  histos.h_eta->Fill(component->eta());
262  histos.h_phi->Fill(component->phi());
263 
264  math::XYZPoint zero(0, 0, 0);
265  math::Error<3>::type zeroCov; // needed for dxyError
266  if (bs) {
267  histos.h_dxy->Fill(component->dxy(*bs));
268  histos.h_exy->Fill(component->dxyError(*bs));
269  } else {
270  histos.h_dxy->Fill(component->dxy(zero));
271  histos.h_exy->Fill(component->dxyError(zero, zeroCov));
272  }
273  if (pv) {
274  histos.h_dz->Fill(component->dz(pv->position()));
275  } else {
276  histos.h_dz->Fill(component->dz(zero));
277  }
278  histos.h_ez->Fill(component->dzError());
279 
280  histos.h_chi2->Fill(component->chi2() / component->ndof());
281 }
ErrorD< N >::type type
Definition: Error.h:32
double pt() const
track transverse momentum
Definition: TrackBase.h:637
double ndof() const
number of degrees of freedom of the fit
Definition: TrackBase.h:590
double dz() const
dz parameter (= dsz/cos(lambda)). This is the track z0 w.r.t (0,0,0) only if the refPoint is close to...
Definition: TrackBase.h:622
double dxyError() const
error on dxy
Definition: TrackBase.h:769
double dzError() const
error on dz
Definition: TrackBase.h:778
def pv(vc)
Definition: MetAnalyzer.py:7
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:649
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:587
histos
Definition: combine.py:4
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:608

◆ fillDecayHistograms()

reco::Vertex const * DiMuonMassBiasMonitor::fillDecayHistograms ( DecayHists const &  histos,
std::vector< const reco::Track *> const &  tracks,
const reco::VertexCollection *const &  pvs,
const edm::EventSetup iSetup 
) const

Definition at line 180 of file DiMuonMassBiasMonitor.cc.

References funct::abs(), TransientTrackBuilder::build(), ChiSquaredProbability(), VertexDistanceXY::distance(), edm::EventSetup::getData(), combine::histos, EgHLTOffHistBins_cfi::mass, mumass2, FSQDQM_cfi::pvs, Measurement1D::significance(), mathSSE::sqrt(), HLT_2022v12_cff::track, tracks, ttbESToken_, Measurement1D::value(), HltBtagValidation_cff::Vertex, KalmanVertexFitter::vertex(), and reco::Vertex::z().

Referenced by analyze().

183  {
184  if (tracks.size() != 2) {
185  edm::LogWarning("DiMuonVertexMonitor") << "There are not enough tracks to construct a vertex!";
186  return nullptr;
187  }
188 
189  const TransientTrackBuilder* theB = &iSetup.getData(ttbESToken_);
190  TransientVertex mumuTransientVtx;
191  std::vector<reco::TransientTrack> tks;
192 
193  for (const auto& track : tracks) {
194  reco::TransientTrack trajectory = theB->build(track);
195  tks.push_back(trajectory);
196  }
197 
198  KalmanVertexFitter kalman(true);
199  mumuTransientVtx = kalman.vertex(tks);
200 
201  auto svtx = reco::Vertex(mumuTransientVtx);
202  if (not svtx.isValid()) {
203  return nullptr;
204  }
205 
206  const auto& mom_t1 = tracks[1]->momentum();
207  const auto& mom_t0 = tracks[0]->momentum();
208  const auto& momentum = mom_t1 + mom_t0;
209 
210  TLorentzVector p4_t0(mom_t0.x(), mom_t0.y(), mom_t0.z(), sqrt(mom_t0.mag2() + mumass2));
211  TLorentzVector p4_t1(mom_t1.x(), mom_t1.y(), mom_t1.z(), sqrt(mom_t1.mag2() + mumass2));
212 
213  const auto& p4 = p4_t0 + p4_t1;
214  float mass = p4.M();
215 
216  auto pvtx = std::min_element(pvs->begin(), pvs->end(), [svtx](reco::Vertex const& pv1, reco::Vertex const& pv2) {
217  return abs(pv1.z() - svtx.z()) < abs(pv2.z() - svtx.z());
218  });
219 
220  if (pvtx == pvs->end()) {
221  return nullptr;
222  }
223 
224  VertexDistanceXY vdistXY;
225  Measurement1D distXY = vdistXY.distance(svtx, *pvtx);
226 
227  auto pvtPos = pvtx->position();
228  const auto& svtPos = svtx.position();
229 
230  math::XYZVector displVect2D(svtPos.x() - pvtPos.x(), svtPos.y() - pvtPos.y(), 0);
231  auto cosAlpha = displVect2D.Dot(momentum) / (displVect2D.Rho() * momentum.rho());
232 
233  auto ct = distXY.value() * cosAlpha * mass / momentum.rho();
234 
235  histos.h_pointing->Fill(cosAlpha);
236 
237  histos.h_mass->Fill(mass);
238 
239  histos.h_pt->Fill(momentum.rho());
240  histos.h_eta->Fill(momentum.eta());
241  histos.h_phi->Fill(momentum.phi());
242 
243  histos.h_ct->Fill(ct);
244 
245  histos.h_displ2D->Fill(distXY.value());
246  histos.h_sign2D->Fill(distXY.significance());
247 
248  if (svtx.chi2() >= 0) {
249  histos.h_vertNormChi2->Fill(svtx.chi2() / svtx.ndof());
250  histos.h_vertProb->Fill(ChiSquaredProbability(svtx.chi2(), svtx.ndof()));
251  }
252 
253  return &*pvtx;
254 }
double z() const
z coordinate
Definition: Vertex.h:133
Measurement1D distance(const GlobalPoint &vtx1Position, const GlobalError &vtx1PositionError, const GlobalPoint &vtx2Position, const GlobalError &vtx2PositionError) const override
static constexpr float mumass2
reco::TransientTrack build(const reco::Track *p) const
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
float ChiSquaredProbability(double chiSquared, double nrDOF)
bool getData(T &iHolder) const
Definition: EventSetup.h:122
auto const & tracks
cannot be loose
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
const edm::ESGetToken< TransientTrackBuilder, TransientTrackRecord > ttbESToken_
double value() const
Definition: Measurement1D.h:25
double significance() const
Definition: Measurement1D.h:29
histos
Definition: combine.py:4
Log< level::Warning, false > LogWarning

◆ fillDescriptions()

void DiMuonMassBiasMonitor::fillDescriptions ( edm::ConfigurationDescriptions descriptions)
static

Definition at line 283 of file DiMuonMassBiasMonitor.cc.

References edm::ParameterSetDescription::add(), edm::ConfigurationDescriptions::addWithDefaultLabel(), submitPVResolutionJobs::desc, HLT_2022v12_cff::InputTag, and AlCaHLTBitMon_QueryRunRegistry::string.

283  {
285  desc.add<edm::InputTag>("muonTracks", edm::InputTag("ALCARECOTkAlDiMuon"));
286  desc.add<edm::InputTag>("vertices", edm::InputTag("offlinePrimaryVertices"));
287  desc.add<edm::InputTag>("beamSpot", edm::InputTag("offlineBeamSpot"));
288  desc.add<std::string>("FolderName", "DiMuonMassBiasMonitor");
289  desc.add<std::string>("decayMotherName", "Z");
290  desc.add<double>("distanceScaleFactor", 0.1);
291 
292  {
293  edm::ParameterSetDescription psDiMuMass;
294  psDiMuMass.add<std::string>("name", "DiMuMass");
295  psDiMuMass.add<std::string>("title", "M(#mu#mu)");
296  psDiMuMass.add<std::string>("yUnits", "[GeV]");
297  psDiMuMass.add<int>("NxBins", 24);
298  psDiMuMass.add<int>("NyBins", 50);
299  psDiMuMass.add<double>("ymin", 65.);
300  psDiMuMass.add<double>("ymax", 115.);
301  desc.add<edm::ParameterSetDescription>("DiMuMassConfig", psDiMuMass);
302  }
303 
304  descriptions.addWithDefaultLabel(desc);
305 }
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ParameterDescriptionBase * add(U const &iLabel, T const &value)

Member Data Documentation

◆ beamSpotToken_

const edm::EDGetTokenT<reco::BeamSpot> DiMuonMassBiasMonitor::beamSpotToken_
private

Definition at line 111 of file DiMuonMassBiasMonitor.h.

Referenced by analyze().

◆ decayMotherName_

const std::string DiMuonMassBiasMonitor::decayMotherName_
private

Definition at line 114 of file DiMuonMassBiasMonitor.h.

Referenced by bookHistograms().

◆ DiMuMassConfiguration_

edm::ParameterSet DiMuonMassBiasMonitor::DiMuMassConfiguration_
private

Definition at line 116 of file DiMuonMassBiasMonitor.h.

Referenced by bookHistograms().

◆ distanceScaleFactor_

const double DiMuonMassBiasMonitor::distanceScaleFactor_
private

Definition at line 115 of file DiMuonMassBiasMonitor.h.

Referenced by bookHistograms().

◆ histosZmm

DecayHists DiMuonMassBiasMonitor::histosZmm
private

Definition at line 121 of file DiMuonMassBiasMonitor.h.

Referenced by analyze(), and bookHistograms().

◆ MEFolderName_

const std::string DiMuonMassBiasMonitor::MEFolderName_
private

Definition at line 113 of file DiMuonMassBiasMonitor.h.

Referenced by bookHistograms().

◆ tracksToken_

const edm::EDGetTokenT<reco::TrackCollection> DiMuonMassBiasMonitor::tracksToken_
private

Definition at line 109 of file DiMuonMassBiasMonitor.h.

Referenced by analyze().

◆ ttbESToken_

const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> DiMuonMassBiasMonitor::ttbESToken_
private

Definition at line 108 of file DiMuonMassBiasMonitor.h.

Referenced by fillDecayHistograms().

◆ vertexToken_

const edm::EDGetTokenT<reco::VertexCollection> DiMuonMassBiasMonitor::vertexToken_
private

Definition at line 110 of file DiMuonMassBiasMonitor.h.

Referenced by analyze().

◆ ZMassPlots

Definition at line 119 of file DiMuonMassBiasMonitor.h.

Referenced by analyze(), and bookHistograms().