Belle II Software development
AnalysisPhase1StudyModule.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 <beast/analysis/modules/AnalysisPhase1StudyModule.h>
10#include <mdst/dataobjects/MCParticle.h>
11#include <generators/SAD/dataobjects/SADMetaHit.h>
12#include <framework/datastore/StoreArray.h>
13#include <framework/logging/Logger.h>
14#include <framework/gearbox/Const.h>
15#include <cmath>
16
17#include <string>
18
19// ROOT
20#include <TMath.h>
21
22using namespace std;
23
24using namespace Belle2;
25
26//using namespace analysis;
27
28//-----------------------------------------------------------------
29// Register the Module
30//-----------------------------------------------------------------
31REG_MODULE(AnalysisPhase1Study);
32
33//-----------------------------------------------------------------
34// Implementation
35//-----------------------------------------------------------------
36
38{
39 // Set module properties
40 setDescription("Study module for BEAST");
41
42}
43
44AnalysisPhase1StudyModule::~AnalysisPhase1StudyModule()
45{
46}
47
48//This module is a histomodule. Any histogram created here will be saved by the HistoManager module
50{
51 h_count = new TH1F("h_count", "", 10, 0., 10.);
52 for (int i = 0; i < 10; i++) {
53 h_prodvtx[i] = new TH2F(TString::Format("h_prodvtx_%d", i), "", 200, -400., 400., 200, 0., 400.);
54 h_decavtx[i] = new TH2F(TString::Format("h_decavtx_%d", i), "", 200, -400., 400., 200, 0., 400.);
55 h_kineticvz[i] = new TH2F(TString::Format("h_kineticvz_%d", i), "", 200, -400., 400., 1000, 0., 10.);
56 h_kineticvz1[i] = new TH2F(TString::Format("h_kineticvz1_%d", i), "", 200, -400., 400., 1000, 0., 10.);
57 h_kineticvz2[i] = new TH2F(TString::Format("h_kineticvz2_%d", i), "", 200, -400., 400., 1000, 0., 10.);
58 h_kineticvz_zoom[i] = new TH2F(TString::Format("h_kineticvz_zoom_%d", i), "", 200, -400., 400., 1000, 0., 0.01);
59 h_Wkineticvz[i] = new TH2F(TString::Format("h_Wkineticvz_%d", i), "", 200, -400., 400., 1000, 0., 10.);
60 h_Wkineticvz_zoom[i] = new TH2F(TString::Format("h_Wkineticvz_zoom_%d", i), "", 200, -400., 400., 1000, 0., 0.01);
61 h_phive[i] = new TH2F(TString::Format("h_phive_%d", i), "", 360, -180., 180., 1000, 0., 10.);
62 h_phivz[i] = new TH2F(TString::Format("h_phivz_%d", i), "", 200, -400., 400., 360, -180., 180.);
63 h_rve[i] = new TH2F(TString::Format("h_rve_%d", i), "", 200, 0., 200., 1000, 0., 10.);
64 h_thetavz[i] = new TH2F(TString::Format("h_thetavz_%d", i), "", 200, -400., 400., 180, 0., 180.);
65 h_kineticvz[i]->Sumw2();
66 h_kineticvz1[i]->Sumw2();
67 h_kineticvz2[i]->Sumw2();
68 h_kineticvz_zoom[i]->Sumw2();
69 h_Wkineticvz[i]->Sumw2();
70 h_Wkineticvz_zoom[i]->Sumw2();
71 h_phive[i]->Sumw2();
72 h_phivz[i]->Sumw2();
73 h_rve[i]->Sumw2();
74 }
75 h_g4_xy = new TH2F("h_g4_xy", "", 100, -5.99, 5.99, 100, -5.99, 5.99);
76 for (int i = 0; i < 2; i++) {
77 h_sad_xy[i] = new TH2F(TString::Format("h_sad_xy_%d", i), "", 100, -5.99, 5.99, 100, -5.99, 5.99);
78 h_sad_sir[i] = new TH1F(TString::Format("h_sad_sir_%d", i), "", 100, -3.99, 3.99);
79 h_z[i] = new TH1F(TString::Format("h_z_%d", i), "", 100, -3.99, 3.99);
80 h_sad_sall[i] = new TH1F(TString::Format("h_sad_sall_%d", i), "", 100, -1499.99, 1499.99);
81 h_sad_sE[i] = new TH2F(TString::Format("h_sad_sE_%d", i), "", 100, -3.99, 3.99, 1000, 0., 10.);
82 h_sad_sraw[i] = new TH1F(TString::Format("h_sad_sraw_%d", i), "", 100, -1499.99, 1499.99);
83 h_sad_E[i] = new TH1F(TString::Format("h_sad_E_%d", i), "", 8000, 0, 7.99);
84 h_sad_s[i] = new TH1F(TString::Format("h_sad_s_%d", i), "", 400, -3.99, 3.99);
85 }
86 h_dpx = new TH1F("h_dpx", "", 1000, -1., 1.);
87 h_dpy = new TH1F("h_dpy", "", 1000, -1., 1.);
88 h_dE = new TH1F("h_dE", "", 1000, -1., 1.);
89 h_px = new TH1F("h_px", "", 10000, -10., 10.);
90 h_py = new TH1F("h_py", "", 10000, -10., 10.);
91 h_pz = new TH1F("h_pz", "", 10000, -10., 10.);
92 h_dx = new TH1F("h_dx", "", 10000, -400., 400.);
93 h_dy = new TH1F("h_dy", "", 10000, -400., 400.);
94 h_dz = new TH1F("h_dz", "", 10000, -400., 400.);
95 h_E = new TH1F("h_E", "", 8000, 0., 7.99);
96 h_P = new TH1F("h_P", "", 8000, 0., 7.99);
97 h_E->Sumw2();
98 h_P->Sumw2();
99}
100
102{
103 B2INFO("AnalysisPhase1StudyModule: Initialize");
104
105 REG_HISTOGRAM
106
107}
108
110{
111}
112
114{
115 double px = 0;
116 double py = 0;
117 double x = 0;
118 double y = 0;
119 double s = 0;
120 //double xraw = 0;
121 //double yraw = 0;
122 //double sraw = 0;
123 double E = 0;
124 double rate = 0;
125 StoreArray<SADMetaHit> sadMetaHits;
126 for (const auto& sadMetaHit : sadMetaHits) {
127 px = sadMetaHit.getpx();
128 py = sadMetaHit.getpy();
129 x = sadMetaHit.getx();
130 y = sadMetaHit.gety();
131 s = sadMetaHit.gets();
132 //xraw = sadMetaHit.getxraw();
133 //yraw = sadMetaHit.getyraw();
134 double sraw = sadMetaHit.getsraw();
135 E = sadMetaHit.getE();
136 rate = sadMetaHit.getrate();
137 h_sad_sir[0]->Fill(s / 100., rate);
138 h_sad_sall[0]->Fill(s / 100., rate);
139 h_sad_sE[0]->Fill(s / 100., E);
140 h_sad_s[0]->Fill(s / 100.);
141 h_sad_E[0]->Fill(E);
142 h_sad_sraw[0]->Fill(sraw / 100., rate);
143 if (-400. < s && s < 400.) {
144 h_sad_xy[1]->Fill(x, y);
145 h_sad_sir[1]->Fill(s / 100., rate);
146 h_sad_sall[1]->Fill(s / 100., rate);
147 h_sad_sraw[1]->Fill(sraw / 100., rate);
148 h_sad_sE[1]->Fill(s / 100., E);
149 h_sad_s[1]->Fill(s / 100.);
150 h_sad_E[1]->Fill(E);
151 }
152 }
153 int counter = 0;
154 StoreArray<MCParticle> mcParticles;
155 for (const auto& mcParticle : mcParticles) { // start loop over all MC particles
156 int PDG = mcParticle.getPDG();
157 //int* Mother = mcParticle.getMother();
158 float Mass = mcParticle.getMass();
159 float Energy = mcParticle.getEnergy();
160 float Kinetic = Energy - Mass;
161 float prodvtx[3];
162 prodvtx[0] = mcParticle.getProductionVertex().X();
163 prodvtx[1] = mcParticle.getProductionVertex().Y();
164 prodvtx[2] = mcParticle.getProductionVertex().Z();
165 float prodr = sqrt(prodvtx[0] * prodvtx[0] + prodvtx[1] * prodvtx[1]);
166 float decavtx[3];
167 decavtx[0] = mcParticle.getDecayVertex().X();
168 decavtx[1] = mcParticle.getDecayVertex().Y();
169 decavtx[2] = mcParticle.getDecayVertex().Z();
170 float decar = sqrt(decavtx[0] * decavtx[0] + decavtx[1] * decavtx[1]);
171 float mom[3];
172 mom[0] = mcParticle.getMomentum().X();
173 mom[1] = mcParticle.getMomentum().Y();
174 mom[2] = mcParticle.getMomentum().Z();
175 float theta = mcParticle.getMomentum().Theta() * TMath::RadToDeg();
176 float phi = mcParticle.getMomentum().Phi() * TMath::RadToDeg();
177 int partID[9] = {0, 0, 0, 0, 0, 0, 0, 0, 0};
178
179 if (PDG == Const::electron.getPDGCode()) partID[0] = 1; //positron
180 else if (PDG == -Const::electron.getPDGCode()) partID[1] = 1; //electron
181 else if (PDG == Const::photon.getPDGCode()) partID[2] = 1; //photon
182 else if (PDG == Const::neutron.getPDGCode()) partID[3] = 1; //neutron
183 else if (PDG == Const::proton.getPDGCode()) partID[4] = 1; //proton
184 else if (PDG == 1000080160) partID[5] = 1; // O
185 else if (PDG == 1000060120) partID[6] = 1; // C
186 else if (PDG == 1000020040) partID[7] = 1; // He
187 else partID[8] = 1;
188 for (int i = 0; i < 8; i++) {
189 if (partID[i] == 1) {
190 h_count->Fill(i);
191 h_prodvtx[i]->Fill(prodvtx[2], prodr);
192 h_decavtx[i]->Fill(decavtx[2], decar);
193 h_kineticvz[i]->Fill(prodvtx[2], Kinetic);
194 h_Wkineticvz[i]->Fill(prodvtx[2], Kinetic, rate);
195 h_kineticvz_zoom[i]->Fill(prodvtx[2], Kinetic);
196 h_Wkineticvz_zoom[i]->Fill(prodvtx[2], Kinetic, rate);
197 h_thetavz[i]->Fill(prodvtx[2], theta);
198 h_phivz[i]->Fill(prodvtx[2], phi);
199 h_phive[i]->Fill(phi, Kinetic);
200 h_rve[i]->Fill(prodr, Kinetic);
201 if (theta == 0) {
202 h_kineticvz1[i]->Fill(prodvtx[2], Kinetic);
203 }
204 if (phi == 0) {
205 h_kineticvz2[i]->Fill(prodvtx[2], Kinetic);
206 }
207 }
208 }
209
210 //h_count->Fill(9);
211 h_prodvtx[8]->Fill(prodvtx[2], prodr);
212 h_decavtx[8]->Fill(decavtx[2], decar);
213 h_kineticvz[8]->Fill(prodvtx[2], Kinetic);
214 h_kineticvz_zoom[8]->Fill(prodvtx[2], Kinetic);
215 h_Wkineticvz[8]->Fill(prodvtx[2], Kinetic, rate);
216 h_Wkineticvz_zoom[8]->Fill(prodvtx[2], Kinetic, rate);
217 h_thetavz[8]->Fill(prodvtx[2], theta);
218 h_phivz[8]->Fill(prodvtx[2], phi);
219 h_phive[8]->Fill(phi, Kinetic);
220 h_rve[8]->Fill(prodr, Kinetic);
221 if (theta == 0) {
222 h_kineticvz1[8]->Fill(prodvtx[2], Kinetic);
223 }
224 if (phi == 0) {
225 h_kineticvz2[8]->Fill(prodvtx[2], Kinetic);
226 }
227 if (counter == 0) {
228 h_count->Fill(9);
229 h_prodvtx[9]->Fill(prodvtx[2], prodr);
230 h_decavtx[9]->Fill(decavtx[2], decar);
231 h_kineticvz[9]->Fill(prodvtx[2], Kinetic);
232 h_kineticvz_zoom[9]->Fill(prodvtx[2], Kinetic);
233 h_Wkineticvz[9]->Fill(prodvtx[2], Kinetic, rate);
234 h_Wkineticvz_zoom[9]->Fill(prodvtx[2], Kinetic, rate);
235 h_thetavz[9]->Fill(prodvtx[2], theta);
236 h_phivz[9]->Fill(prodvtx[2], phi);
237
238 h_dx->Fill(x - prodvtx[0]);
239 h_dy->Fill(-y - prodvtx[1]);
240 h_dz->Fill(-s - prodvtx[2]);
241
242 h_z[0]->Fill(prodvtx[2] / 100.);
243 h_z[1]->Fill(prodvtx[2] / 100., rate);
244
245 h_dpx->Fill(px - mom[0]);
246 h_dpy->Fill(-py - mom[1]);
247 h_dE->Fill(E - Energy);
248
249 h_g4_xy->Fill(prodvtx[0], prodvtx[1]);
250
251 h_px->Fill(mom[0]);
252 h_py->Fill(mom[1]);
253 h_pz->Fill(mom[2]);
254 h_E->Fill(Kinetic);
255 h_P->Fill(sqrt(mom[0] * mom[0] + mom[1] * mom[1] + mom[2] * mom[2]));
256 }
257 counter++;
258 }
259
260}
261
263{
264
265
266}
267
269{
270}
271
272
R E
internal precision of FFTW codelets
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
virtual void endRun() override
End-of-run action.
virtual void terminate() override
Termination action.
virtual void beginRun() override
Called when entering a new run.
AnalysisPhase1StudyModule()
Constructor: Sets the description, the properties and the parameters of the module.
virtual void defineHisto() override
Defines the histograms.
static const ParticleType neutron
neutron particle
Definition: Const.h:675
static const ChargedStable proton
proton particle
Definition: Const.h:663
static const ParticleType photon
photon particle
Definition: Const.h:673
static const ChargedStable electron
electron particle
Definition: Const.h:659
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
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:649
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.
STL namespace.