Belle II Software development
ARICHNtupleModule.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// Own header.
10#include <arich/modules/arichNtuple/ARICHNtupleModule.h>
11
12// Hit classes
13#include <framework/dataobjects/EventMetaData.h>
14#include <arich/dataobjects/ARICHLikelihood.h>
15#include <arich/dataobjects/ARICHAeroHit.h>
16#include <arich/dataobjects/ARICHTrack.h>
17#include <arich/dataobjects/ARICHPhoton.h>
18#include <arich/dataobjects/ARICHHit.h>
19#include <arich/dataobjects/ARICHInfo.h>
20
21// framework - DataStore
22#include <framework/datastore/StoreObjPtr.h>
23
24#include <mdst/dataobjects/HitPatternCDC.h>
25
26#ifdef ALIGNMENT_USING_BHABHA
27#include <mdst/dataobjects/ECLCluster.h>
28#endif
29
30// framework aux
31#include <framework/gearbox/Const.h>
32#include <framework/logging/Logger.h>
33
34// ROOT
35#include <Math/Vector3D.h>
36#include <vector>
37
38namespace Belle2 {
44 //-----------------------------------------------------------------
46 //-----------------------------------------------------------------
47
48 REG_MODULE(ARICHNtuple);
49
50 //-----------------------------------------------------------------
51 // Implementation
52 //-----------------------------------------------------------------
53
55 m_file(0),
56 m_tree(0)
57 {
58 // set module description (e.g. insert text)
59 setDescription("The module saves variables needed for performance analysis, such as postion and momentum of the hit, likelihoods for hypotheses and number of photons.");
60
61 // Add parameters
62 addParam("outputFile", m_outputFile, "ROOT output file name", std::string("ARICHNtuple.root"));
63 }
64
66 {
67 }
68
70 {
71
72 if (m_file) delete m_file;
73 m_file = new TFile(m_outputFile.c_str(), "RECREATE");
74 if (m_file->IsZombie()) {
75 B2FATAL("Couldn't open file '" << m_outputFile << "' for writing!");
76 return;
77 }
78
79 m_tree = new TTree("arich", "ARICH validation ntuple");
80
81 m_tree->Branch("evt", &m_arich.evt, "evt/I");
82 m_tree->Branch("run", &m_arich.run, "run/I");
83 m_tree->Branch("exp", &m_arich.exp, "exp/I");
84
85 m_tree->Branch("charge", &m_arich.charge, "charge/S");
86 m_tree->Branch("pValue", &m_arich.pValue, "pValue/F");
87 m_tree->Branch("d0", &m_arich.z0, "d0/F");
88 m_tree->Branch("z0", &m_arich.d0, "z0/F");
89
90#ifdef ALIGNMENT_USING_BHABHA
91 m_tree->Branch("eop", &m_arich.eop, "eop/F");
92 m_tree->Branch("e9e21", &m_arich.e9e21, "e9e21/F");
93 m_tree->Branch("etot", &m_arich.etot, "etot/F");
94#endif
95
96 m_tree->Branch("PDG", &m_arich.PDG, "PDG/I");
97 m_tree->Branch("motherPDG", &m_arich.motherPDG, "motherPDG/I");
98 m_tree->Branch("primary", &m_arich.primary, "primary/S");
99 m_tree->Branch("seen", &m_arich.seen, "seen/S");
100 m_tree->Branch("scatter", &m_arich.scatter, "scatter/I");
101
102 m_tree->Branch("rhoProd", &m_arich.rhoProd, "rhoProd/F");
103 m_tree->Branch("zProd", &m_arich.zProd, "zProd/F");
104 m_tree->Branch("phiProd", &m_arich.phiProd, "phiProd/F");
105 m_tree->Branch("rhoDec", &m_arich.rhoDec, "rhoDec/F");
106 m_tree->Branch("zDec", &m_arich.zDec, "zDec/F");
107 m_tree->Branch("phiDec", &m_arich.phiDec, "phiDec/F");
108 m_tree->Branch("status", &m_arich.status, "status/I");
109
110 m_tree->Branch("detPhot", &m_arich.detPhot, "detPhot/I");
111 m_tree->Branch("numBkg", &m_arich.numBkg, "e/F:mu:pi:K:p:d");
112 m_tree->Branch("expPhot", &m_arich.expPhot, "e/F:mu:pi:K:p:d");
113 m_tree->Branch("logL", &m_arich.logL, "e/F:mu:pi:K:p:d");
114
115 m_tree->Branch("recHit", &m_arich.recHit, "PDG/I:x/F:y:z:p:theta:phi");
116 m_tree->Branch("mcHit", &m_arich.mcHit, "PDG/I:x/F:y:z:p:theta:phi");
117 m_tree->Branch("winHit", &m_arich.winHit, "x/F:y");
118 m_tree->Branch("nrec", &m_arich.nRec, "nRec/I");
119 m_tree->Branch("nCDC", &m_arich.nCDC, "nCDC/I");
120 m_tree->Branch("inAcceptance", &m_arich.inAcc, "inAcc/O");
121 m_tree->Branch("photons", "std::vector<Belle2::ARICHPhoton>", &m_arich.photons);
122
123 // required input
124 m_arichTracks.isRequired();
125 m_arichLikelihoods.isRequired();
126
127 // optional input
128 m_tracks.isOptional();
130 m_arichAeroHits.isOptional();
131 m_arichInfo.isOptional();
132
133 StoreArray<ARICHHit> arichHits;
134 arichHits.isRequired();
135
136 }
137
139 {
140
141 StoreObjPtr<EventMetaData> evtMetaData;
142 StoreArray<ARICHHit> arichHits;
143
144 if (!m_arichTracks.isValid()) return;
145
146#ifdef ALIGNMENT_USING_BHABHA
147 double E_tot = 0.;
148 StoreArray<ECLCluster> eclClusters;
149 for (const auto& trk : m_tracks) {
150 const ECLCluster* eclTrack = trk.getRelated<ECLCluster>();
151 if (eclTrack) {
152 if (eclTrack->getEnergy() > 0.1) E_tot += eclTrack->getEnergy();
153 }
154 }
155 for (const auto& cluster : eclClusters) {
156 if (!cluster.isNeutral()) continue;
157 if (cluster.getHypothesisId() != 5 && cluster.getHypothesisId() != 1)
158 continue;
159 if (cluster.getEnergy() > 0.1) E_tot += cluster.getEnergy();
160 }
161#endif
162
163 for (const auto& arichTrack : m_arichTracks) {
164
165 const ExtHit* extHit = arichTrack.getRelated<ExtHit>();
166
167 const Track* mdstTrack = NULL;
168 if (extHit) mdstTrack = extHit->getRelated<Track>();
169
170 const ARICHAeroHit* aeroHit = arichTrack.getRelated<ARICHAeroHit>();
171
172 const ARICHLikelihood* lkh = NULL;
173 if (mdstTrack) lkh = mdstTrack->getRelated<ARICHLikelihood>();
174 else lkh = arichTrack.getRelated<ARICHLikelihood>();
175
176 if (!lkh) continue;
177
178 m_arich.clear();
179
180 m_arich.evt = evtMetaData->getEvent();
181 m_arich.run = evtMetaData->getRun();
182 m_arich.exp = evtMetaData->getExperiment();
183
184 // set hapd window hit if available
185 if (arichTrack.hitsWindow()) {
186 m_arich.winHit[0] = arichTrack.windowHitPosition().X();
187 m_arich.winHit[1] = arichTrack.windowHitPosition().Y();
188 }
189
190 if (lkh->getFlag() == 1) m_arich.inAcc = 1;
191 m_arich.logL.e = lkh->getLogL(Const::electron);
192 m_arich.logL.mu = lkh->getLogL(Const::muon);
193 m_arich.logL.pi = lkh->getLogL(Const::pion);
194 m_arich.logL.K = lkh->getLogL(Const::kaon);
195 m_arich.logL.p = lkh->getLogL(Const::proton);
196 m_arich.logL.d = lkh->getLogL(Const::deuteron);
197
198 m_arich.expPhot.e = lkh->getExpPhot(Const::electron);
199 m_arich.expPhot.mu = lkh->getExpPhot(Const::muon);
200 m_arich.expPhot.pi = lkh->getExpPhot(Const::pion);
201 m_arich.expPhot.K = lkh->getExpPhot(Const::kaon);
202 m_arich.expPhot.p = lkh->getExpPhot(Const::proton);
203 m_arich.expPhot.d = lkh->getExpPhot(Const::deuteron);
204
205 m_arich.detPhot = lkh->getDetPhot();
206
207#ifdef ALIGNMENT_USING_BHABHA
208 m_arich.etot = E_tot;
209#endif
210 m_arich.status = 1;
211
212 m_arich.trgtype = 0;
213 if (m_arichInfo.getEntries() > 0) {
214 auto arichInfo = *m_arichInfo[0];
215 m_arich.trgtype = arichInfo.gettrgtype();
216 }
217
218 const MCParticle* particle = 0;
219
220 if (mdstTrack) {
221 const TrackFitResult* fitResult = mdstTrack->getTrackFitResultWithClosestMass(Const::pion);
222 if (fitResult) {
223 m_arich.pValue = fitResult->getPValue();
224 ROOT::Math::XYZVector trkPos = fitResult->getPosition();
225 m_arich.charge = fitResult->getChargeSign();
226 m_arich.z0 = trkPos.Z();
227 m_arich.d0 = trkPos.Rho();
228 m_arich.nCDC = fitResult->getHitPatternCDC().getNHits();
229#ifdef ALIGNMENT_USING_BHABHA
230 ROOT::Math::XYZVector trkMom = fitResult->getMomentum();
231 const ECLCluster* eclTrack = mdstTrack->getRelated<ECLCluster>();
232 if (eclTrack) {
233 m_arich.eop = eclTrack->getEnergy() / trkMom.R();
234 m_arich.e9e21 = eclTrack->getE9oE21();
235 }
236#endif
237 }
238 m_arich.status += 10;
239 int fromARICHMCP = 0;
240 particle = mdstTrack->getRelated<MCParticle>();
241 if (!particle) { particle = mdstTrack->getRelated<MCParticle>("arichMCParticles"); fromARICHMCP = 1;}
242 if (particle) {
243 m_arich.PDG = particle->getPDG();
244 m_arich.primary = particle->getStatus(MCParticle::c_PrimaryParticle);
245 m_arich.seen = particle->hasSeenInDetector(Const::ARICH);
246 ROOT::Math::XYZVector prodVertex = particle->getProductionVertex();
247 m_arich.rhoProd = prodVertex.Rho();
248 m_arich.zProd = prodVertex.Z();
249 m_arich.phiProd = prodVertex.Phi();
250 ROOT::Math::XYZVector decVertex = particle->getDecayVertex();
251 m_arich.rhoDec = decVertex.Rho();
252 m_arich.zDec = decVertex.Z();
253 m_arich.phiDec = decVertex.Phi();
254
255 if (!fromARICHMCP) {
256 MCParticle* mother = particle->getMother();
257 if (mother) m_arich.motherPDG = mother->getPDG();
258 std::vector<Belle2::MCParticle*> daughs = particle->getDaughters();
259 for (const auto daugh : daughs) {
260 if (daugh == NULL) continue;
261 if (daugh->getPDG() == particle->getPDG()) m_arich.scatter = 1;
262 }
263 }
264 }
265 }
266
267
268 // get reconstructed photons associated with track
269 const std::vector<ARICHPhoton>& photons = arichTrack.getPhotons();
270 m_arich.nRec = photons.size();
271 for (auto it = photons.begin(); it != photons.end(); ++it) {
272 ARICHPhoton iph = *it;
273 if (iph.getHitID() < arichHits.getEntries()) {
274 auto ihit = *arichHits[iph.getHitID()];
275 int module = ihit.getModule();
276 int channel = ihit.getChannel();
277 int chid = 1 << 16 | module << 8 | channel;
278 iph.setHitID(chid);
279 }
280 m_arich.photons.push_back(iph);
281 }
282
283 ROOT::Math::XYZVector recPos = arichTrack.getPosition();
284 m_arich.recHit.x = recPos.X();
285 m_arich.recHit.y = recPos.Y();
286 m_arich.recHit.z = recPos.Z();
287
288 ROOT::Math::XYZVector recMom = arichTrack.getDirection() * arichTrack.getMomentum();
289 m_arich.recHit.p = recMom.R();
290 m_arich.recHit.theta = recMom.Theta();
291 m_arich.recHit.phi = recMom.Phi();
292
293 if (aeroHit) {
294 ROOT::Math::XYZVector truePos = aeroHit->getPosition();
295 m_arich.mcHit.x = truePos.X();
296 m_arich.mcHit.y = truePos.Y();
297 m_arich.mcHit.z = truePos.Z();
298
299 ROOT::Math::XYZVector trueMom = aeroHit->getMomentum();
300 m_arich.mcHit.p = trueMom.R();
301 m_arich.mcHit.theta = trueMom.Theta();
302 m_arich.mcHit.phi = trueMom.Phi();
303
304 m_arich.mcHit.PDG = aeroHit->getPDG();
305 m_arich.status += 1000;
306
307 if (!particle) {
308 particle = aeroHit->getRelated<MCParticle>();
309 if (particle) {
310 m_arich.PDG = particle->getPDG();
311 MCParticle* mother = particle->getMother();
312 if (mother) m_arich.motherPDG = mother->getPDG();
313 m_arich.primary = particle->getStatus(MCParticle::c_PrimaryParticle);
314 m_arich.seen = particle->hasSeenInDetector(Const::ARICH);
315 ROOT::Math::XYZVector prodVertex = particle->getProductionVertex();
316 m_arich.rhoProd = prodVertex.Rho();
317 m_arich.zProd = prodVertex.Z();
318 m_arich.phiProd = prodVertex.Phi();
319 ROOT::Math::XYZVector decVertex = particle->getDecayVertex();
320 m_arich.rhoDec = decVertex.Rho();
321 m_arich.zDec = decVertex.Z();
322 m_arich.phiDec = decVertex.Phi();
323
324 std::vector<Belle2::MCParticle*> daughs = particle->getDaughters();
325 for (const auto daugh : daughs) {
326 if (daugh->getPDG() == particle->getPDG()) m_arich.scatter = 1;
327 }
328 }
329 }
330 }
331 m_tree->Fill();
332 }
333 }
334
336 {
337 m_file->cd();
338 m_tree->Write();
339 m_file->Close();
340 }
341
343} // end Belle2 namespace
344
Datastore class that holds information on track parameters at the entrance in aerogel.
Definition: ARICHAeroHit.h:27
This is a class to store ARICH likelihoods in the datastore.
TTree * m_tree
pointer to output tree
StoreArray< ARICHAeroHit > m_arichAeroHits
Optional input array of ARICHAeroHits.
std::string m_outputFile
output root file
StoreArray< ARICHTrack > m_arichTracks
Required array of input ARICHTracks.
ARICH::ARICHTree m_arich
ntuple structure
StoreArray< MCParticle > m_arichMCPs
Optional input array of MCParticles.
StoreArray< Track > m_tracks
Optional input array of Tracks.
StoreArray< ARICHLikelihood > m_arichLikelihoods
Required array of input ARICHLikelihoods.
TFile * m_file
pointer to output root file
StoreArray< ARICHInfo > m_arichInfo
Optional input array of ARICHInfo.
Struct for ARICH reconstructed photon (hit related to track) information.
Definition: ARICHPhoton.h:25
void setHitID(int id)
Set ID of corresponding ARICHHit.
Definition: ARICHPhoton.h:48
int getHitID()
Get ID of corresponding ARICHHit.
Definition: ARICHPhoton.h:230
static const ChargedStable muon
muon particle
Definition: Const.h:660
static const ChargedStable pion
charged pion particle
Definition: Const.h:661
static const ChargedStable proton
proton particle
Definition: Const.h:663
static const ChargedStable kaon
charged kaon particle
Definition: Const.h:662
static const ChargedStable electron
electron particle
Definition: Const.h:659
static const ChargedStable deuteron
deuteron particle
Definition: Const.h:664
ECL cluster data.
Definition: ECLCluster.h:27
double getE9oE21() const
Return E9/E21 (shower shape variable).
Definition: ECLCluster.h:280
double getEnergy(EHypothesisBit hypothesis) const
Return Energy (GeV).
Definition: ECLCluster.cc:23
Store one Ext hit as a ROOT object.
Definition: ExtHit.h:31
unsigned short getNHits() const
Get the total Number of CDC hits in the fit.
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:47
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
Values of the result of a track fit with a given particle hypothesis.
short getChargeSign() const
Return track charge (1 or -1).
double getPValue() const
Getter for Chi2 Probability of the track fit.
ROOT::Math::XYZVector getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
ROOT::Math::XYZVector getPosition() const
Getter for vector of position at closest approach of track in r/phi projection.
HitPatternCDC getHitPatternCDC() const
Getter for the hit pattern in the CDC;.
Class that bundles various TrackFitResults.
Definition: Track.h:25
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
virtual void terminate() override
Termination action.
virtual ~ARICHNtupleModule()
Destructor.
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#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.
Float_t zDec
Decay vertex (cylindrical coordinate z) of MCParticle.
Int_t detPhot
Number of detected photons.
Int_t nCDC
Number of track CDC hits.
Float_t winHit[2]
Track hit on hapd window (x,y coordinates).
Float_t rhoDec
Decay vertex (cylindrical coordinate r) of MCParticle.
ParticlesArray logL
Log likelihoods.
Int_t scatter
1 if particle scattered (i.e.
Int_t nRec
Number of reconstructed photons.
Short_t primary
Is a primary particle (from related MCParticle).
Int_t evt
Event number.
Short_t seen
Is seen in ARICH (from related MCParticle).
TrackHit recHit
Extrapolated Track hit.
Int_t motherPDG
PDG code of related mother MCParticle.
Int_t exp
Experiment number.
ParticlesArray numBkg
Number of expected background photons.
ParticlesArray expPhot
Number of expected photons (signal + bkg).
std::vector< Belle2::ARICHPhoton > photons
Vector of reconstructed photons.
Float_t phiProd
Production vertex (cylindrical coordinate phi) of MCParticle.
Int_t trgtype
Event trigger type.
Float_t rhoProd
Production vertex (cylindrical coordinate r) of MCParticle.
void clear()
Clear the structure: set elements to zero.
Float_t phiDec
Decay vertex (cylindrical coordinate phi) of MCParticle.
Bool_t inAcc
Track in detector acceptance, i.e.
TrackHit mcHit
Related MC particle hit.
Int_t PDG
PDG code of related MCParticle.
Int_t status
Track status (add proper description).
Float_t pValue
p-value of Track fit.
Float_t zProd
Production vertex (cylindrical coordinate z) of MCParticle.
Float_t p
momentum magnitude
Float_t z
impact point, z component
Float_t theta
momentum polar angle
Float_t phi
momentum azimuthal angle
Float_t x
impact point, x component
Float_t y
impact point, y component