Belle II Software development
TrackAna.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 <daq/dqm/modules/TrackAna.h>
10
11#include <framework/dataobjects/EventMetaData.h>
12#include <framework/datastore/StoreObjPtr.h>
13#include <framework/datastore/StoreArray.h>
14#include <mdst/dataobjects/Track.h>
15#include <mdst/dataobjects/TrackFitResult.h>
16
17using namespace std;
18using namespace Belle2;
19
20//-----------------------------------------------------------------
21// Register the Module
22//-----------------------------------------------------------------
23REG_MODULE(TrackAna);
24
25//-----------------------------------------------------------------
26// Implementation
27//-----------------------------------------------------------------
28
30{
31 //Set module properties
32 setDescription("The simplest physics analysis");
34
35 //Parameter definition
36 B2INFO("TrackAna: Constructor done.");
37}
38
39
40TrackAnaModule::~TrackAnaModule()
41{
42}
43
45{
46 REG_HISTOGRAM
47}
48
50{
51 h_multi = new TH1F("Multi", "Particle Multiplicity", 50, 0.0, 50.0);
52 h_p[0] = new TH1F("Px", "Particle Momentum X", 100, -5.0, 5.0);
53 h_p[1] = new TH1F("Py", "Particle Momentum Y", 100, -5.0, 5.0);
54 h_p[2] = new TH1F("Pz", "Particle Momentum Z", 100, -5.0, 5.0);
55 h_p[3] = new TH1F("E", "Particle Energy", 100, 0.0, 10.0);
56}
57
59{
60 B2INFO("TrackAna: started to measure elapsed time.");
61}
62
63
65{
66 // Get Event Info
68
69 // Get List of Tracks
70 StoreArray<Track> trklist;
71
72 int ntrk = trklist.getEntries();
73 h_multi->Fill((float)ntrk);
74 // printf ( "***** Generator Event List : Event %d ; multi = %d *****\n", evtno, npart );
75 for (int i = 0; i < ntrk; i++) {
76 Track* trk = trklist[i];
78 ROOT::Math::PxPyPzEVector p4 = fit->get4Momentum();
79 h_p[0]->Fill(p4.Px());
80 h_p[1]->Fill(p4.Py());
81 h_p[2]->Fill(p4.Pz());
82 h_p[3]->Fill(p4.E());
83 }
84}
85
87{
88}
89
90
92{
93 B2INFO("TrackAna: terminate called");
94}
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
void initialize() override
Module functions to be called from main process.
Definition: TrackAna.cc:44
void event() override
Function to process event record.
Definition: TrackAna.cc:64
void endRun() override
Function to process end_run record.
Definition: TrackAna.cc:86
void terminate() override
Function to terminate module.
Definition: TrackAna.cc:91
void beginRun() override
Module functions to be called from event process.
Definition: TrackAna.cc:58
TrackAnaModule()
Constructor / Destructor.
Definition: TrackAna.cc:29
void defineHisto() override
Module funcions to define histograms.
Definition: TrackAna.cc:49
Values of the result of a track fit with a given particle hypothesis.
Class that bundles various TrackFitResults.
Definition: Track.h:25
const TrackFitResult * getTrackFitResult(const Const::ChargedStable &chargedStable) const
Default Access to TrackFitResults.
Definition: Track.cc:30
#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.
STL namespace.