Belle II Software
release-05-01-25
BeamSpotCollectorModule.cc
1
/**************************************************************************
2
* BASF2 (Belle Analysis Framework 2) *
3
* Copyright(C) 2017 - Belle II Collaboration *
4
* *
5
* Author: The Belle II Collaboration *
6
* Contributors: Radek Zlebcik
7
* *
8
* *
9
* This software is provided "as is" without any warranty. *
10
**************************************************************************/
11
#include <tracking/modules/BeamSpotCollector/BeamSpotCollectorModule.h>
12
13
#include <analysis/dataobjects/ParticleList.h>
14
#include <analysis/utility/ReferenceFrame.h>
15
#include <mdst/dataobjects/TrackFitResult.h>
16
17
using namespace
Belle2
;
18
using namespace
std;
19
20
#include <iostream>
21
22
//-----------------------------------------------------------------
23
// Register the Module
24
//-----------------------------------------------------------------
25
REG_MODULE
(BeamSpotCollector)
26
27
//-----------------------------------------------------------------
28
// Implementation
29
//-----------------------------------------------------------------
30
31
BeamSpotCollectorModule
::
BeamSpotCollectorModule
() :
CalibrationCollectorModule
(),
32
m_evt(-99), m_exp(-99), m_run(-99),
33
m_time(-99),
34
m_mu0_d0(-99), m_mu0_z0(-99), m_mu0_phi0(-99), m_mu0_tanlambda(-99), m_mu0_omega(-99),
35
m_mu1_d0(-99), m_mu1_z0(-99), m_mu1_phi0(-99), m_mu1_tanlambda(-99), m_mu1_omega(-99)
36
{
37
//Set module properties
38
39
setDescription(
"Collect data for BeamSpot calibration algorithm using the position of mu+mu- events"
);
40
setPropertyFlags(c_ParallelProcessingCertified);
41
42
addParam(
"Y4SPListName"
, m_Y4SPListName,
"Name of the Y4S particle list"
, std::string(
"Upsilon(4S):IPDQM"
));
43
}
44
45
void
BeamSpotCollectorModule::prepare
()
46
{
47
B2INFO(
"Init of the trees"
);
48
std::string objectName =
"events"
;
49
//Data object creation --------------------------------------------------
50
TTree* tree =
new
TTree(objectName.c_str(),
""
);
51
52
tree->Branch<
int
>(
"event"
, &m_evt);
53
tree->Branch<
int
>(
"exp"
, &m_exp);
54
tree->Branch<
int
>(
"run"
, &m_run);
55
tree->Branch<
double
>(
"time"
, &m_time);
56
57
tree->Branch<
double
>(
"mu0_d0"
, &m_mu0_d0);
58
tree->Branch<
double
>(
"mu0_z0"
, &m_mu0_z0);
59
tree->Branch<
double
>(
"mu0_phi0"
, &m_mu0_phi0);
60
tree->Branch<
double
>(
"mu0_tanlambda"
, &m_mu0_tanlambda);
61
62
tree->Branch<
double
>(
"mu1_d0"
, &m_mu1_d0);
63
tree->Branch<
double
>(
"mu1_z0"
, &m_mu1_z0);
64
tree->Branch<
double
>(
"mu1_phi0"
, &m_mu1_phi0);
65
tree->Branch<
double
>(
"mu1_tanlambda"
, &m_mu1_tanlambda);
66
67
68
// We register the objects so that our framework knows about them.
69
// Don't try and hold onto the pointers or fill these objects directly
70
// Use the getObjectPtr functions to access collector objects
71
registerObject<TTree>(objectName, tree);
72
}
73
74
75
void
BeamSpotCollectorModule::collect
()
76
{
77
m_evt = m_emd->getEvent();
78
m_run = m_emd->getRun();
79
m_exp = m_emd->getExperiment();
80
m_time = m_emd->getTime() / 1e9 / 3600.;
//from ns to hours
81
82
83
StoreObjPtr<ParticleList>
Y4SParticles(m_Y4SPListName);
84
85
86
if
(!Y4SParticles.
isValid
() || abs(Y4SParticles->getPDGCode()) != 300553)
87
return
;
88
89
if
(Y4SParticles->getListSize() != 1)
90
return
;
91
92
93
std::vector<int> indxes = Y4SParticles->getParticle(0)->getDaughterIndices();
94
if
(indxes.size() != 2)
return
;
95
96
const
Particle
* part0 = Y4SParticles->getParticle(0)->getDaughter(0);
97
const
Particle
* part1 = Y4SParticles->getParticle(0)->getDaughter(1);
98
const
TrackFitResult
* tr0 = part0->
getTrackFitResult
();
99
const
TrackFitResult
* tr1 = part1->
getTrackFitResult
();
100
101
m_mu0_d0 = tr0->
getD0
();
102
m_mu0_z0 = tr0->
getZ0
();
103
m_mu0_phi0 = tr0->
getPhi0
();
104
m_mu0_tanlambda = tr0->
getTanLambda
();
105
m_mu0_omega = tr0->
getOmega
();
106
107
108
m_mu1_d0 = tr1->
getD0
();
109
m_mu1_z0 = tr1->
getZ0
();
110
m_mu1_phi0 = tr1->
getPhi0
();
111
m_mu1_tanlambda = tr1->
getTanLambda
();
112
m_mu1_omega = tr1->
getOmega
();
113
114
115
getObjectPtr<TTree>(
"events"
)->Fill();
116
117
}
Belle2::CalibrationCollectorModule
Calibration collector module base class.
Definition:
CalibrationCollectorModule.h:44
Belle2::BeamSpotCollectorModule
This collects the track parameters of the mu+mu- events for calibration of the BeamSpot using CAF and...
Definition:
BeamSpotCollectorModule.h:33
Belle2::Particle::getTrackFitResult
const TrackFitResult * getTrackFitResult() const
Returns the pointer to the TrackFitResult that was used to create this Particle (ParticleType == c_Tr...
Definition:
Particle.cc:779
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition:
Module.h:652
Belle2::TrackFitResult::getOmega
double getOmega() const
Getter for omega.
Definition:
TrackFitResult.h:194
Belle2::BeamSpotCollectorModule::collect
void collect() override final
Event processor The filling of the tree.
Definition:
BeamSpotCollectorModule.cc:75
Belle2::TrackFitResult::getTanLambda
double getTanLambda() const
Getter for tanLambda.
Definition:
TrackFitResult.h:206
Belle2::TrackFitResult
Values of the result of a track fit with a given particle hypothesis.
Definition:
TrackFitResult.h:59
Belle2::TrackFitResult::getZ0
double getZ0() const
Getter for z0.
Definition:
TrackFitResult.h:200
Belle2
Abstract base class for different kinds of events.
Definition:
MillepedeAlgorithm.h:19
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition:
ParticleList.h:33
Belle2::TrackFitResult::getPhi0
double getPhi0() const
Getter for phi0.
Definition:
TrackFitResult.h:184
Belle2::BeamSpotCollectorModule::prepare
void prepare() override final
Initialize the module.
Definition:
BeamSpotCollectorModule.cc:45
Belle2::Particle
Class to store reconstructed particles.
Definition:
Particle.h:77
Belle2::TrackFitResult::getD0
double getD0() const
Getter for d0.
Definition:
TrackFitResult.h:178
Belle2::StoreObjPtr::isValid
bool isValid() const
Check whether the object was created.
Definition:
StoreObjPtr.h:120
tracking
modules
BeamSpotCollector
src
BeamSpotCollectorModule.cc
Generated on Fri Nov 5 2021 04:00:32 for Belle II Software by
1.8.17