Belle II Software  release-06-00-14
PhysicsObjectsMiraBelleDst2Module.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #include <dqm/modules/PhysicsObjectsMiraBelle/PhysicsObjectsMiraBelleDst2Module.h>
10 #include <analysis/dataobjects/ParticleList.h>
11 #include <analysis/variables/ContinuumSuppressionVariables.h>
12 #include <analysis/variables/TrackVariables.h>
13 #include <analysis/utility/PCmsLabTransform.h>
14 #include <analysis/variables/Variables.h>
15 #include <framework/datastore/StoreObjPtr.h>
16 #include <framework/datastore/StoreArray.h>
17 #include <mdst/dataobjects/Track.h>
18 #include <mdst/dataobjects/TrackFitResult.h>
19 #include <mdst/dataobjects/KLMCluster.h>
20 #include <mdst/dataobjects/HitPatternCDC.h>
21 #include <mdst/dataobjects/HitPatternVXD.h>
22 #include <mdst/dataobjects/EventLevelTrackingInfo.h>
23 #include <mdst/dataobjects/PIDLikelihood.h>
24 #include <top/variables/TOPDigitVariables.h>
25 #include <arich/modules/arichDQM/ARICHDQMModule.h>
26 #include <arich/dataobjects/ARICHLikelihood.h>
27 #include <mdst/dataobjects/SoftwareTriggerResult.h>
28 #include <TDirectory.h>
29 #include <map>
30 
31 using namespace Belle2;
32 
33 REG_MODULE(PhysicsObjectsMiraBelleDst2)
34 
35 PhysicsObjectsMiraBelleDst2Module::PhysicsObjectsMiraBelleDst2Module() : HistoModule()
36 {
37  setDescription("Monitor Physics Objects Quality");
38  setPropertyFlags(c_ParallelProcessingCertified);
39 
40  addParam("TriggerIdentifier", m_triggerIdentifier,
41  "Trigger identifier string used to select events for the histograms", std::string("software_trigger_cut&skim&accept_hadron"));
42  addParam("DstListName", m_dstListName, "Name of the D*+ particle list", std::string("D*+:physMiraBelle"));
43 }
44 
45 void PhysicsObjectsMiraBelleDst2Module::defineHisto()
46 {
47  TDirectory* oldDir = gDirectory;
48  oldDir->mkdir("PhysicsObjectsMiraBelleDst2")->cd();
49  // Mass distributions
50  m_h_D0_pi0_InvM = new TH1F("hist_D0_pi0_InvM", "#pi^{0} mass of D^{*}#rightarrowK#pi#pi^{0};D0_pi0_InvM;", 50, 0.09, 0.17);
51 
52  oldDir->cd();
53 }
54 
55 
57 {
58  REG_HISTOGRAM
59 
61  result.isOptional();
62 }
63 
65 {
66  m_h_D0_pi0_InvM->Reset();
67 }
68 
70 {
71 
73  if (!result.isValid()) {
74  B2WARNING("SoftwareTriggerResult object not available but needed to select events for the histograms.");
75  return;
76  }
77 
78  const std::map<std::string, int>& results = result->getResults();
79  if (results.find(m_triggerIdentifier) == results.end()) {
80  B2WARNING("PhysicsObjectsMiraBelleDst2: Can't find trigger identifier: " << m_triggerIdentifier);
81  return;
82  }
83 
84  // apply software trigger
85  const bool accepted = (result->getResult(m_triggerIdentifier) == SoftwareTriggerCutResult::c_accept);
86  if (accepted == false) return;
87 
88  // get D* candidates
90 
91  for (unsigned int i = 0; i < dstParticles->getListSize(); i++) {
92  const Particle* dst = dstParticles->getParticle(i);
93  const Particle* d0 = dst->getDaughter(0);
94  float dst_mass = dst->getMass();
95  float d0_mass = d0->getMass();
96  float delta_m = dst_mass - d0_mass;
97  // Mass of Pi0
98  if (1.83 < d0_mass && d0_mass < 1.89 && 0.142 < delta_m && delta_m < 0.148) {
99  m_h_D0_pi0_InvM->Fill(Belle2::Variable::particleInvariantMassFromDaughters(d0->getDaughter(2)));
100  }
101  }
102 }
103 
104 // void PhysicsObjectsMiraBelleDst2Module::getPIDInformationMiraBelle(Particle* part, float pid[3][7])
105 // {
106 // PIDLikelihood* I_like_coffee = part->getPIDLikelihood();
107 // pid[0][0] = I_like_coffee->getProbability();
108 // }
109 
111 {
112 }
113 
115 {
116 }
117 
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
Class to store reconstructed particles.
Definition: Particle.h:74
float getMass() const
Returns invariant mass (= nominal for FS particles)
Definition: Particle.h:445
const Particle * getDaughter(unsigned i) const
Returns a pointer to the i-th daughter particle.
Definition: Particle.cc:639
std::string m_dstListName
Name of the mu+ particle list.
void initialize() override
Function for dynamic initialization of module.
void event() override
Function to process event record.
void endRun() override
Function to process end_run record.
void terminate() override
Function to terminate module.
std::string m_triggerIdentifier
Trigger identifier string used to select events for the histograms.
void beginRun() override
Function to process begin_run record.
TH1F * m_h_D0_pi0_InvM
Pi0 invariant mass for D0->Kpipi0.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
@ c_accept
Accept this event.
Abstract base class for different kinds of events.