Belle II Software  release-05-02-19
BeamSpotMonitorModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2019 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Giulia Casarosa, Gaetano De Marino *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <tracking/modules/ipMonitor/BeamSpotMonitorModule.h>
12 #include <framework/datastore/StoreObjPtr.h>
13 #include <framework/dataobjects/EventMetaData.h>
14 
15 #include <TVector3.h>
16 #include <TMatrixDSym.h>
17 
18 using namespace Belle2;
19 using namespace std;
20 
21 //-----------------------------------------------------------------
22 // Register the Module
23 //-----------------------------------------------------------------
24 REG_MODULE(BeamSpotMonitor)
25 
26 //-----------------------------------------------------------------
27 // Implementation
28 //-----------------------------------------------------------------
29 
31 {
32  // Set module properties
33  setDescription("Module for the monitoring of the BeamSpot position and size");
34 
35  // Parameter definitions
36  addParam("outputFileName", m_rootFileName, "Name of output root file.", std::string("BeamSpotMonitor.root"));
37 }
38 
40 {
41 
42  TDirectory* olddir = gDirectory;
43  m_rootFilePtr = new TFile(m_rootFileName.c_str(), "RECREATE");
44 
45  //tree initialization
46  m_tree = new TTree("bspt", "RECREATE");
47  b_exp = m_tree->Branch("exp", &m_exp, "exp/i");
48  b_run = m_tree->Branch("run", &m_run, "run/i");
49  b_x = m_tree->Branch("x", &m_x, "x/d");
50  b_y = m_tree->Branch("y", &m_y, "y/d");
51  b_z = m_tree->Branch("z", &m_z, "z/d");
52  b_xErr = m_tree->Branch("xErr", &m_xErr, "xErr/d");
53  b_yErr = m_tree->Branch("yErr", &m_yErr, "yErr/d");
54  b_zErr = m_tree->Branch("zErr", &m_zErr, "zErr/d");
55  b_xSize = m_tree->Branch("xSize", &m_xSize, "xSize/d");
56  b_ySize = m_tree->Branch("ySize", &m_ySize, "ySize/d");
57  b_zSize = m_tree->Branch("zSize", &m_zSize, "zSize/d");
58 
59  olddir->cd();
60 
61 }
62 
64 {
65  if (! m_BeamSpotDB.isValid()) {
66  B2WARNING("No valid BeamSpot for the requested IoV");
67  } else {
68  m_BeamSpot = *m_BeamSpotDB;
69  }
70 
71 }
72 
74 {
75 
77  m_exp = meta->getExperiment();
78  m_run = meta->getRun();
79  B2DEBUG(25, "monitoring beam spot for experiment = " << m_exp << ", run = " << m_run);
80 
81  if (! m_BeamSpotDB.isValid())
82  return;
83 
84  //retrieve vertex position
85  m_x = m_BeamSpot.getIPPosition().X();
86  m_y = m_BeamSpot.getIPPosition().Y();
87  m_z = m_BeamSpot.getIPPosition().Z();
88 
89  //retrieve vertex position error
90  TMatrixDSym posErr = m_BeamSpot.getIPPositionCovMatrix();
91  m_xErr = sqrt(posErr[0][0]);
92  m_yErr = sqrt(posErr[1][1]);
93  m_zErr = sqrt(posErr[2][2]);
94 
95  //retrieve beam spot size
96  TMatrixDSym size = m_BeamSpot.getSizeCovMatrix();
97  m_xSize = sqrt(size[0][0]);
98  m_ySize = sqrt(size[1][1]);
99  m_zSize = sqrt(size[2][2]);
100 
101  m_tree->Fill();
102 
103 }
104 
106 {
107  TDirectory* olddir = gDirectory;
108  if (m_rootFilePtr != nullptr) {
109 
110  m_rootFilePtr->cd();
111 
112  //write the tree
113  m_tree->Write();
114 
115  m_rootFilePtr->Close();
116 
117  olddir->cd();
118  }
119 }
120 
Belle2::BeamSpotMonitorModule
Module for the monitoring of the BeamSpot position and size.
Definition: BeamSpotMonitorModule.h:39
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::BeamSpotMonitorModule::terminate
virtual void terminate() override
print the payloads uniqueID and write tree to the rootfile
Definition: BeamSpotMonitorModule.cc:105
Belle2::BeamSpotMonitorModule::initialize
virtual void initialize() override
initialize the TTree
Definition: BeamSpotMonitorModule.cc:39
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::BeamSpotMonitorModule::beginRun
virtual void beginRun() override
check BeamSpot payload validity
Definition: BeamSpotMonitorModule.cc:63
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::BeamSpotMonitorModule::event
virtual void event() override
fill trees
Definition: BeamSpotMonitorModule.cc:73