Belle II Software  release-06-00-14
MasterClassModule.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 <masterclass/modules/MasterClassModule.h>
10 #include <mdst/dataobjects/PIDLikelihood.h>
11 
12 #include "TDirectory.h"
13 
14 using namespace Belle2;
15 
16 //-----------------------------------------------------------------
17 // Register the Module
18 //-----------------------------------------------------------------
19 REG_MODULE(MasterClass)
20 
21 //-----------------------------------------------------------------
22 // Implementation
23 //-----------------------------------------------------------------
24 
26 {
27  // Set module properties
28  setDescription(R"DOC("Module to write out data in a format for Belle II masterclasses)DOC");
29  addParam("outputFileName", m_filename, "File name of root ntuple output.", std::string("masterclass.root"));
30 }
31 
33 {
34  m_tracks.isRequired();
35  m_clusters.isRequired();
36 
37  m_event = new BEvent;
38  m_file = TFile::Open(m_filename.c_str(), "RECREATE");
39  m_tree = new TTree("T", "Event Tree");
40  m_tree->Branch("BEvent", &m_event);
41 }
42 
44 {
45  m_event->Clear();
47 
48  for (auto& track : m_tracks) {
49  auto pid = track.getRelated<PIDLikelihood>();
50  const double priors[] = {0.05, 0.05, 0.65, 0.24, 0.01, 0};
51  auto type = pid->getMostLikely(priors);
52  auto trackFit = track.getTrackFitResultWithClosestMass(type);
53  auto p = trackFit->getMomentum();
54  double m = type.getMass();
55  SIMPLEPID id = ALL;
56  switch (type.getPDGCode()) {
57  case 11: id = ELECTRON; break;
58  case 13: id = MUON; break;
59  case 211: id = PION; break;
60  case 321: id = KAON; break;
61  case 2212: id = PROTON; break;
62  default: id = ALL;
63  }
64  m_event->AddTrack(p.x(), p.y(), p.z(), sqrt(m * m + p.Mag2()), trackFit->getChargeSign(), id);
65  }
66 
67  for (auto& cluster : m_clusters) {
68  if (!cluster.hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) continue;
69  double E = cluster.getEnergy(ECLCluster::EHypothesisBit::c_nPhotons);
70  if (E < 0.1) continue;
71  TVector3 p = cluster.getClusterPosition();
72  p.SetMag(E);
73  m_event->AddTrack(p.x(), p.y(), p.z(), E, 0, PHOTON);
74  }
75 
76  m_tree->Fill();
77 }
78 
80 {
81  m_tree->Write();
82  m_file->Close();
83  delete m_file;
84 }
The Class for Masterclass event parameters.
Definition: BEvent.h:19
virtual void Clear(Option_t *="")
Clear the array of particles.
Definition: BEvent.cc:50
void AddTrack(float px, float py, float pz, float e, float charge, SIMPLEPID pid)
Add the track to the event.
Definition: BEvent.cc:33
void EventNo(int evtno)
Set the current event number.
Definition: BEvent.cc:23
@ c_nPhotons
CR is split into n photons (N1)
Module to write out data in a format for Belle II masterclasses.
TTree * m_tree
output tree
virtual void terminate() override
Close ntuple.
BEvent * m_event
output event object
StoreArray< ECLCluster > m_clusters
Cluster objects.
virtual void event() override
Write out particles.
StoreArray< Track > m_tracks
Track objects.
virtual void initialize() override
Register input and output data.
TFile * m_file
root ntuple file
std::string m_filename
output file name
Base class for Modules.
Definition: Module.h:72
Class to collect log likelihoods from TOP, ARICH, dEdx, ECL and KLM aimed for output to mdst includes...
Definition: PIDLikelihood.h:26
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.