10 #include <arich/modules/arichNtuple/ARICHNtupleModule.h>
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>
22 #include <framework/datastore/StoreObjPtr.h>
24 #include <mdst/dataobjects/HitPatternCDC.h>
26 #ifdef ALIGNMENT_USING_BHABHA
27 #include <mdst/dataobjects/ECLCluster.h>
31 #include <framework/gearbox/Const.h>
32 #include <framework/logging/Logger.h>
61 setDescription(
"The module saves variables needed for performance analysis, such as postion and momentum of the hit, likelihoods for hypotheses and number of photons.");
64 addParam(
"outputFile", m_outputFile,
"ROOT output file name",
string(
"ARICHNtuple.root"));
67 ARICHNtupleModule::~ARICHNtupleModule()
71 void ARICHNtupleModule::initialize()
74 if (m_file)
delete m_file;
75 m_file =
new TFile(m_outputFile.c_str(),
"RECREATE");
76 if (m_file->IsZombie()) {
77 B2FATAL(
"Couldn't open file '" << m_outputFile <<
"' for writing!");
81 m_tree =
new TTree(
"arich",
"ARICH validation ntuple");
83 m_tree->Branch(
"evt", &m_arich.evt,
"evt/I");
84 m_tree->Branch(
"run", &m_arich.run,
"run/I");
85 m_tree->Branch(
"exp", &m_arich.exp,
"exp/I");
87 m_tree->Branch(
"charge", &m_arich.charge,
"charge/S");
88 m_tree->Branch(
"pValue", &m_arich.pValue,
"pValue/F");
89 m_tree->Branch(
"d0", &m_arich.z0,
"pValue/F");
90 m_tree->Branch(
"z0", &m_arich.d0,
"pValue/F");
92 #ifdef ALIGNMENT_USING_BHABHA
93 m_tree->Branch(
"eop", &m_arich.eop,
"eop/F");
94 m_tree->Branch(
"e9e21", &m_arich.e9e21,
"e9e21/F");
95 m_tree->Branch(
"etot", &m_arich.etot,
"etot/F");
98 m_tree->Branch(
"PDG", &m_arich.PDG,
"PDG/I");
99 m_tree->Branch(
"motherPDG", &m_arich.motherPDG,
"motherPDG/I");
100 m_tree->Branch(
"primary", &m_arich.primary,
"primary/S");
101 m_tree->Branch(
"seen", &m_arich.seen,
"seen/S");
102 m_tree->Branch(
"scatter", &m_arich.scatter,
"scatter/I");
104 m_tree->Branch(
"rhoProd", &m_arich.rhoProd,
"rhoProd/F");
105 m_tree->Branch(
"zProd", &m_arich.zProd,
"zProd/F");
106 m_tree->Branch(
"phiProd", &m_arich.phiProd,
"phiProd/F");
107 m_tree->Branch(
"rhoDec", &m_arich.rhoDec,
"rhoDec/F");
108 m_tree->Branch(
"zDec", &m_arich.zDec,
"zDec/F");
109 m_tree->Branch(
"phiDec", &m_arich.phiDec,
"phiDec/F");
110 m_tree->Branch(
"status", &m_arich.status,
"status/I");
112 m_tree->Branch(
"detPhot", &m_arich.detPhot,
"detPhot/I");
113 m_tree->Branch(
"numBkg", &m_arich.numBkg,
"e/F:mu:pi:K:p:d");
114 m_tree->Branch(
"expPhot", &m_arich.expPhot,
"e/F:mu:pi:K:p:d");
115 m_tree->Branch(
"logL", &m_arich.logL,
"e/F:mu:pi:K:p:d");
117 m_tree->Branch(
"recHit", &m_arich.recHit,
"PDG/I:x/F:y:z:p:theta:phi");
118 m_tree->Branch(
"mcHit", &m_arich.mcHit,
"PDG/I:x/F:y:z:p:theta:phi");
119 m_tree->Branch(
"winHit", &m_arich.winHit,
"x/F:y");
120 m_tree->Branch(
"nrec", &m_arich.nRec,
"nRec/I");
121 m_tree->Branch(
"nCDC", &m_arich.nCDC,
"nCDC/I");
122 m_tree->Branch(
"inAcceptance", &m_arich.inAcc,
"inAcc/O");
123 m_tree->Branch(
"photons",
"std::vector<Belle2::ARICHPhoton>", &m_arich.photons);
126 m_arichTracks.isRequired();
127 m_arichLikelihoods.isRequired();
130 m_tracks.isOptional();
131 m_arichMCPs.isOptional();
132 m_arichAeroHits.isOptional();
133 m_arichInfo.isOptional();
136 arichHits.isRequired();
140 void ARICHNtupleModule::beginRun()
144 void ARICHNtupleModule::event()
150 if (!m_arichTracks.isValid())
return;
152 #ifdef ALIGNMENT_USING_BHABHA
155 for (
const auto& trk : m_tracks) {
161 for (
const auto& cluster : eclClusters) {
162 if (!cluster.isNeutral())
continue;
163 if (cluster.getHypothesisId() != 5 && cluster.getHypothesisId() != 1)
165 if (cluster.getEnergy() > 0.1) E_tot += cluster.getEnergy();
169 for (
const auto& arichTrack : m_arichTracks) {
173 const Track* mdstTrack = NULL;
186 m_arich.evt = evtMetaData->getEvent();
187 m_arich.run = evtMetaData->getRun();
188 m_arich.exp = evtMetaData->getExperiment();
191 if (arichTrack.hitsWindow()) {
192 TVector2 winHit = arichTrack.windowHitPosition();
193 m_arich.winHit[0] = winHit.X();
194 m_arich.winHit[1] = winHit.Y();
197 if (lkh->getFlag() == 1) m_arich.inAcc = 1;
198 m_arich.logL.e = lkh->getLogL(Const::electron);
199 m_arich.logL.mu = lkh->getLogL(Const::muon);
200 m_arich.logL.pi = lkh->getLogL(Const::pion);
201 m_arich.logL.K = lkh->getLogL(Const::kaon);
202 m_arich.logL.p = lkh->getLogL(Const::proton);
203 m_arich.logL.d = lkh->getLogL(Const::deuteron);
205 m_arich.expPhot.e = lkh->getExpPhot(Const::electron);
206 m_arich.expPhot.mu = lkh->getExpPhot(Const::muon);
207 m_arich.expPhot.pi = lkh->getExpPhot(Const::pion);
208 m_arich.expPhot.K = lkh->getExpPhot(Const::kaon);
209 m_arich.expPhot.p = lkh->getExpPhot(Const::proton);
210 m_arich.expPhot.d = lkh->getExpPhot(Const::deuteron);
212 m_arich.detPhot = lkh->getDetPhot();
214 #ifdef ALIGNMENT_USING_BHABHA
215 m_arich.etot = E_tot;
220 if (m_arichInfo.getEntries() > 0) {
221 auto arichInfo = *m_arichInfo[0];
222 m_arich.trgtype = arichInfo.gettrgtype();
228 const TrackFitResult* fitResult = mdstTrack->getTrackFitResultWithClosestMass(Const::pion);
233 m_arich.z0 = trkPos.Z();
234 m_arich.d0 = (trkPos.XYvector()).Mod();
236 #ifdef ALIGNMENT_USING_BHABHA
240 m_arich.eop = eclTrack->
getEnergy() / trkMom.Mag();
245 m_arich.status += 10;
246 int fromARICHMCP = 0;
247 particle = mdstTrack->getRelated<
MCParticle>();
248 if (!particle) { particle = mdstTrack->getRelated<
MCParticle>(
"arichMCParticles"); fromARICHMCP = 1;}
250 m_arich.PDG = particle->getPDG();
251 m_arich.primary = particle->getStatus(MCParticle::c_PrimaryParticle);
252 m_arich.seen = particle->hasSeenInDetector(Const::ARICH);
253 TVector3 prodVertex = particle->getProductionVertex();
254 m_arich.rhoProd = prodVertex.Perp();
255 m_arich.zProd = prodVertex.Z();
256 m_arich.phiProd = prodVertex.Phi();
257 TVector3 decVertex = particle->getDecayVertex();
258 m_arich.rhoDec = decVertex.Perp();
259 m_arich.zDec = decVertex.Z();
260 m_arich.phiDec = decVertex.Phi();
264 if (mother) m_arich.motherPDG = mother->getPDG();
265 std::vector<Belle2::MCParticle*> daughs = particle->getDaughters();
266 for (
const auto daugh : daughs) {
267 if (daugh == NULL)
continue;
268 if (daugh->getPDG() == particle->getPDG()) m_arich.scatter = 1;
276 const std::vector<ARICHPhoton>& photons = arichTrack.getPhotons();
277 m_arich.nRec = photons.size();
279 for (
auto it = photons.begin(); it != photons.end(); ++it) {
281 if (iph.
getHitID() < arichHits.getEntries()) {
282 auto ihit = *arichHits[iph.
getHitID()];
283 int module = ihit.getModule();
284 int channel = ihit.getChannel();
285 int chid = 1 << 16 | module << 8 | channel;
288 m_arich.photons.push_back(iph);
292 TVector3 recPos = arichTrack.getPosition();
293 m_arich.recHit.x = recPos.X();
294 m_arich.recHit.y = recPos.Y();
295 m_arich.recHit.z = recPos.Z();
297 TVector3 recMom = arichTrack.getDirection() * arichTrack.getMomentum();
298 m_arich.recHit.p = recMom.Mag();
299 m_arich.recHit.theta = recMom.Theta();
300 m_arich.recHit.phi = recMom.Phi();
303 TVector3 truePos = aeroHit->getPosition();
304 m_arich.mcHit.x = truePos.X();
305 m_arich.mcHit.y = truePos.Y();
306 m_arich.mcHit.z = truePos.Z();
308 TVector3 trueMom = aeroHit->getMomentum();
309 m_arich.mcHit.p = trueMom.Mag();
310 m_arich.mcHit.theta = trueMom.Theta();
311 m_arich.mcHit.phi = trueMom.Phi();
313 m_arich.mcHit.PDG = aeroHit->getPDG();
314 m_arich.status += 1000;
319 m_arich.PDG = particle->getPDG();
321 if (mother) m_arich.motherPDG = mother->getPDG();
322 m_arich.primary = particle->getStatus(MCParticle::c_PrimaryParticle);
323 m_arich.seen = particle->hasSeenInDetector(Const::ARICH);
324 TVector3 prodVertex = particle->getProductionVertex();
325 m_arich.rhoProd = prodVertex.Perp();
326 m_arich.zProd = prodVertex.Z();
327 m_arich.phiProd = prodVertex.Phi();
328 TVector3 decVertex = particle->getDecayVertex();
329 m_arich.rhoDec = decVertex.Perp();
330 m_arich.zDec = decVertex.Z();
331 m_arich.phiDec = decVertex.Phi();
333 std::vector<Belle2::MCParticle*> daughs = particle->getDaughters();
334 for (
const auto daugh : daughs) {
335 if (daugh->getPDG() == particle->getPDG()) m_arich.scatter = 1;
346 void ARICHNtupleModule::endRun()
350 void ARICHNtupleModule::terminate()
Datastore class that holds information on track parameters at the entrance in aerogel.
This is a class to store ARICH likelihoods in the datastore.
ARICH reconstruction efficiency test module.
Struct for ARICH reconstructed photon (hit related to track) information.
void setHitID(int id)
Set ID of corresponding ARICHHit.
int getHitID()
Get ID of corresponding ARICHHit.
double getE9oE21() const
Return E9/E21 (shower shape variable).
double getEnergy(EHypothesisBit hypothesis) const
Return Energy (GeV).
Store one Ext hit as a ROOT object.
unsigned short getNHits() const
Get the total Number of CDC hits in the fit.
A Class to store the Monte Carlo particle information.
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
Accessor to arrays stored in the data store.
Type-safe access to single objects in the data store.
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.
TVector3 getMomentum() const
Getter for vector of momentum at closest approach of track in r/phi projection.
TVector3 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.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.