Belle II Software development
HitLevelInfoWriter.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 <cdc/modules/HitLevelInfoWriter/HitLevelInfoWriter.h>
10
11#include <analysis/dataobjects/ParticleList.h>
12#include <analysis/dataobjects/Particle.h>
13#include <analysis/utility/ReferenceFrame.h>
14#include <mdst/dataobjects/HitPatternCDC.h>
15#include <mdst/dataobjects/PIDLikelihood.h>
16#include <mdst/dbobjects/BeamSpot.h>
17
18#include <TMath.h>
19
20using namespace Belle2;
21using namespace std;
22
23REG_MODULE(HitLevelInfoWriter);
24
26{
27
28 setDescription("Extract dE/dx information for calibration development.");
29
30 addParam("outputBaseName", m_strOutputBaseName, "Suffix for output file name", string("HLInfo.root"));
31 addParam("particleLists", m_strParticleList, "Vector of ParticleLists to save", vector<string>());
32 addParam("enableHitLevel", m_isHitLevel, "True or False for Hit level variables", false);
33 addParam("m_isExtraVar", m_isExtraVar, "True or False for extra track/hit level variables", false);
34 addParam("nodeadwire", m_isDeadwire, "True or False for deadwire hit variables", false);
35 addParam("relativeCorrections", m_isRelative, "If true, apply corrections relative to those used in production", false);
36 addParam("corrections", m_isCorrection, "If true, apply corrections", true);
37
38}
39
41
43{
44
45 B2INFO("Creating a ROOT file for the hit level information...");
46
47 // required inputs
48 m_dedxTracks.isRequired();
49 m_tracks.isRequired();
50 m_trackFitResults.isRequired();
51 m_eclClusters.isOptional();
52 m_klmClusters.isOptional();
53
54 // build a map to relate input strings to the right particle type
55 map<string, string> pdgMap = {{"pi+", "Const::pion.getPDGCode()"}, {"K+", "Const::kaon.getPDGCode()"}, {"mu+", "Const::muon.getPDGCode()"}, {"e+", "Const::electron.getPDGCode()"}, {"p+", "Const::proton.getPDGCode()"}, {"deuteron", "Const::deuteron.getPDGCode()"}};
56
57 // if no particle lists are given, write out all tracks
59
60 // create a new output file for each particle list specified
61 for (unsigned int i = 0; i < m_strParticleList.size(); i++) {
62 // strip the name of the particle lists to make this work
63 string pdg = pdgMap[m_strParticleList[i].substr(0, m_strParticleList[i].find(":"))];
64 string filename = string(m_strOutputBaseName + "_PID" + pdg + ".root");
65 bookOutput(filename);
66 }
67}
68
70{
71
72 int nParticleList = m_strParticleList.size();
73
74 // **************************************************
75 // LOOP OVER dE/dx measurements for all tracks if
76 // no particle list is specified.
77 // **************************************************
78 if (nParticleList == 0) {
79 for (int idedx = 0; idedx < m_dedxTracks.getEntries(); idedx++) {
80
81 CDCDedxTrack* dedxTrack = m_dedxTracks[idedx];
82 if (!dedxTrack) {
83 B2WARNING("No dedx related track...");
84 continue;
85 }
86
87 const Track* track = dedxTrack->getRelatedFrom<Track>();
88 if (!track) {
89 B2WARNING("No related track...");
90 continue;
91 }
92
93 const TrackFitResult* fitResult = track->getTrackFitResultWithClosestMass(Const::pion);
94 if (!fitResult) {
95 B2WARNING("No related fit for this track...");
96 continue;
97 }
98
99 if (dedxTrack->size() == 0 || dedxTrack->size() > 200) continue;
100 if (dedxTrack->getCosTheta() < TMath::Cos(150.0 * TMath::DegToRad()))continue; //-0.866
101 if (dedxTrack->getCosTheta() > TMath::Cos(17.0 * TMath::DegToRad())) continue; //0.95
102
103 // fill the event meta data
104 StoreObjPtr<EventMetaData> evtMetaData;
105 m_expID = evtMetaData->getExperiment();
106 m_runID = evtMetaData->getRun();
107 m_eventID = evtMetaData->getEvent();
108
109 m_injring = dedxTrack->getInjectionRing();
110 m_injtime = dedxTrack->getInjectionTime();
111
112 // fill the E/P
113 const ECLCluster* eclCluster = track->getRelated<ECLCluster>();
114 if (eclCluster and eclCluster->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) {
115 m_eop = (eclCluster->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons)) / (fitResult->getMomentum().R());
117 if (m_isExtraVar) {
118 m_e1_9 = eclCluster->getE1oE9();
119 m_e9_21 = eclCluster->getE9oE21();
120 m_eclsnHits = eclCluster->getNumberOfCrystals();
121 }
122 // fill the muon depth
123 const KLMCluster* klmCluster = eclCluster->getRelated<KLMCluster>();
124 if (klmCluster) m_klmLayers = klmCluster->getLayers();
125 }
126
127 // fill the TTree with the Track information
128 fillTrack(fitResult);
129
130 // fill the TTree with the CDCDedxTrack information
131 fillDedx(dedxTrack);
132
133 // fill the TTree
134 m_tree[0]->Fill();
135 }
136 }
137
138 // **************************************************
139 //
140 // LOOP OVER particles in the given particle lists
141 //
142 // **************************************************
143 for (int iList = 0; iList < nParticleList; iList++) {
144 // make sure the list exists and is not empty
145 StoreObjPtr<ParticleList> particlelist(m_strParticleList[iList]);
146 if (!particlelist or particlelist->getListSize(true) == 0) {
147 //B2WARNING("ParticleList " << m_strParticleList[iList] << " not found or empty, skipping");
148 continue;
149 }
150
151 // loop over the particles in the list and follow the links to the
152 // dE/dx information (Particle -> PIDLikelihood -> Track -> CDCDedxTrack)
153 for (unsigned int iPart = 0; iPart < particlelist->getListSize(true); iPart++) {
154 Particle* part = particlelist->getParticle(iPart, true);
155 if (!part) {
156 B2WARNING("No particles...");
157 continue;
158 }
160 if (!pid) {
161 B2WARNING("No related PID likelihood...");
162 continue;
163 }
164 Track* track = pid->getRelatedFrom<Track>();
165 if (!track) {
166 B2WARNING("No related track...");
167 continue;
168 }
169 CDCDedxTrack* dedxTrack = track->getRelatedTo<CDCDedxTrack>();
170 if (!dedxTrack) {
171 B2WARNING("No related CDCDedxTrack...");
172 continue;
173 }
174 string ptype = m_strParticleList[iList].substr(0, m_strParticleList[iList].find(":"));
175 const TrackFitResult* fitResult = track->getTrackFitResult(Const::pion);
176 if (ptype != "pi+") {
177 if (ptype == "K+") fitResult = track->getTrackFitResultWithClosestMass(Const::kaon);
178 else if (ptype == "p+") fitResult = track->getTrackFitResultWithClosestMass(Const::proton);
179 else if (ptype == "deuteron") fitResult = track->getTrackFitResultWithClosestMass(Const::deuteron);
180 else if (ptype == "mu+") fitResult = track->getTrackFitResultWithClosestMass(Const::muon);
181 else if (ptype == "e+") fitResult = track->getTrackFitResultWithClosestMass(Const::electron);
182 }
183 if (!fitResult) {
184 B2WARNING("No related fit for this track...");
185 continue;
186 }
187
188 if (dedxTrack->size() == 0 || dedxTrack->size() > 200) continue;
189 //if out CDC (dont add as we dont correct via correctionmodules)
190 if (dedxTrack->getCosTheta() < TMath::Cos(150.0 * TMath::DegToRad()))continue; //-0.866
191 if (dedxTrack->getCosTheta() > TMath::Cos(17.0 * TMath::DegToRad())) continue; //0.95
192
193 // fill the event meta data
194 StoreObjPtr<EventMetaData> evtMetaData;
195 m_expID = evtMetaData->getExperiment();
196 m_runID = evtMetaData->getRun();
197 m_eventID = evtMetaData->getEvent();
198
199 m_injring = dedxTrack->getInjectionRing();
200 m_injtime = dedxTrack->getInjectionTime();
201
202 // fill the E/P
203 const ECLCluster* eclCluster = track->getRelated<ECLCluster>();
204 if (eclCluster and eclCluster->hasHypothesis(ECLCluster::EHypothesisBit::c_nPhotons)) {
205 m_eop = (eclCluster->getEnergy(ECLCluster::EHypothesisBit::c_nPhotons)) / (fitResult->getMomentum().R());
207 if (m_isExtraVar) {
208 m_e1_9 = eclCluster->getE1oE9();
209 m_e9_21 = eclCluster->getE9oE21();
210 m_eclsnHits = eclCluster->getNumberOfCrystals();
211 }
212 // fill the muon depth
213 const KLMCluster* klmCluster = eclCluster->getRelated<KLMCluster>();
214 if (klmCluster) m_klmLayers = klmCluster->getLayers();
215 }
216
217 // fill the TTree with the Track information
218 fillTrack(fitResult);
219
220 // fill the TTree with the CDCDedxTrack information
221 fillDedx(dedxTrack);
222
223 // fill the TTree
224 m_tree[iList]->Fill();
225 }
226 }
227}
228
229//---------------------------------------------------------------------------------------
231{
232
233 for (unsigned int i = 0; i < m_file.size(); i++) {
234 B2INFO("Done writing out the hit level information...\t" << m_tree[i]->GetEntries() << " tracks");
235
236 // write the ttree to a root file
237 m_file[i]->cd();
238 m_tree[i]->Write();
239 m_file[i]->Close();
240 }
241}
242
243//---------------------------------------------------------------------------------------
245{
246 ROOT::Math::XYZVector trackMom = fitResult->getMomentum();
247 m_p = trackMom.R();
248 m_pt = trackMom.Rho();
249 m_phi = trackMom.Phi();
250
251 m_theta = trackMom.Theta() * 180. / TMath::Pi(); //in degree
252 if (m_theta > 17. && m_theta < 150.)m_inCDC = 1;
253 else m_inCDC = 0;
254
255 if (fitResult->getChargeSign() < 0) {
256 m_p *= -1;
257 m_pt *= -1;
258 }
259
260 ROOT::Math::XYZVector trackPos = fitResult->getPosition();
261 m_vx0 = trackPos.X();
262 m_vy0 = trackPos.Y();
263 m_vz0 = trackPos.Z();
264
265 m_d0 = fitResult->getD0();
266 m_z0 = fitResult->getZ0();
267 m_chi2 = fitResult->getPValue();
268 m_tanlambda = fitResult->getTanLambda();
269 m_phi0 = fitResult->getPhi0();
270 m_nCDChits = fitResult->getHitPatternCDC().getNHits();
271
272 static DBObjPtr<BeamSpot> beamSpotDB;
273 const auto& frame = ReferenceFrame::GetCurrent();
274 UncertainHelix helix = fitResult->getUncertainHelix();
275 helix.passiveMoveBy(ROOT::Math::XYZVector(beamSpotDB->getIPPosition()));
276 m_dr = frame.getVertex(ROOT::Math::XYZVector(helix.getPerigee())).Rho();
277 m_dphi = frame.getVertex(ROOT::Math::XYZVector(helix.getPerigee())).Phi();
278 m_dz = frame.getVertex(ROOT::Math::XYZVector(helix.getPerigee())).Z();
279
280}
281
282//---------------------------------------------------------------------------------------
284{
285 // clear the containers first
286 clearEntries();
287
288 m_trackID = dedxTrack->trackID();
289 m_length = dedxTrack->getLength();
290 m_charge = dedxTrack->getCharge();
291 m_cosTheta = dedxTrack->getCosTheta();
292 m_PDG = dedxTrack->getPDG();
293
294 m_pCDC = dedxTrack->getMomentum();
295 if (m_charge < 0) {
296 m_pCDC *= -1;
297 }
298
299 h_nhits = dedxTrack->size();
300
301 // Get the calibration constants
302 m_scale = m_DBScaleFactor->getScaleFactor();
303 m_runGain = m_DBRunGain->getRunGain();
304 m_timeGain = m_DBInjectTime->getCorrection("mean", m_injring, m_injtime);
305 m_timeReso = m_DBInjectTime->getCorrection("reso", m_injring, m_injtime);
306
309 } else {
310 m_cosEdgeCor = 1.0;
311 }
312 m_hadronpars = m_DBHadronCor->getHadronPars();
313
314 //variable to save layer variables
315 map<int, vector<double>>l_var;
316 double cdcChi[Const::ChargedStable::c_SetSize];
317
318 //Modify the hit level dedx
319 if (m_isCorrection) {
320 recalculateDedx(dedxTrack, l_var, cdcChi);
321 m_chie = cdcChi[0];
322 m_chimu = cdcChi[1];
323 m_chipi = cdcChi[2];
324 m_chik = cdcChi[3];
325 m_chip = cdcChi[4];
326 m_chid = cdcChi[5];
327 } else {
328 m_mean = dedxTrack->getDedxMean();
329 m_trunc = dedxTrack->getDedx();
330 m_truncNoSat = dedxTrack->getDedxNoSat();
331 m_error = dedxTrack->getDedxError();
332 m_chie = dedxTrack->getChi(0);
333 m_chimu = dedxTrack->getChi(1);
334 m_chipi = dedxTrack->getChi(2);
335 m_chik = dedxTrack->getChi(3);
336 m_chip = dedxTrack->getChi(4);
337 m_chid = dedxTrack->getChi(5);
338 }
339
340 m_pmeane = dedxTrack->getPmean(0);
341 m_pmeanmu = dedxTrack->getPmean(1);
342 m_pmeanpi = dedxTrack->getPmean(2);
343 m_pmeank = dedxTrack->getPmean(3);
344 m_pmeanp = dedxTrack->getPmean(4);
345 m_pmeand = dedxTrack->getPmean(5);
346
347 m_prese = dedxTrack->getPreso(0);
348 m_presmu = dedxTrack->getPreso(1);
349 m_prespi = dedxTrack->getPreso(2);
350 m_presk = dedxTrack->getPreso(3);
351 m_presp = dedxTrack->getPreso(4);
352 m_presd = dedxTrack->getPreso(5);
353
354 // Get the vector of dE/dx values for all layers
355 double lout = 0, lin = 0, increment = 0;
356 int lastlayer = 0;
357
358 if (m_isCorrection) {
359 l_nhits = l_var[0].size();
360 const int lEdgeTrunc = int(l_nhits * 0.05 + 0.51);
361 const int hEdgeTrunc = int(l_nhits * (1 - 0.25) + 0.51);
362 l_nhitsused = hEdgeTrunc - lEdgeTrunc;
363 } else {
364 l_nhits = dedxTrack->getNLayerHits();
365 l_nhitsused = dedxTrack->getNLayerHitsUsed();
366 }
367
368 for (int il = 0; il < l_nhits; ++il) {
369 if (m_isCorrection) {
370 l_nhitscombined[il] = l_var[0][il];
371 l_wirelongesthit[il] = l_var[1][il];
372 l_layer[il] = l_var[2][il];
373 l_path[il] = l_var[3][il];
374 l_dedx[il] = l_var[4][il];
375 } else {
376 l_nhitscombined[il] = dedxTrack->getNHitsCombined(il);
377 l_wirelongesthit[il] = dedxTrack->getWireLongestHit(il);
378 l_layer[il] = dedxTrack->getLayer(il);
379 l_path[il] = dedxTrack->getLayerPath(il);
380 l_dedx[il] = dedxTrack->getLayerDedx(il);
381 }
382 l_cosCor[il] = m_DBCosineCor->getMean(l_layer[il], m_cosTheta);
383
384 if (l_layer[il] > lastlayer) lout++;
385 else if (l_layer[il] < lastlayer) lin++;
386 else continue;
387
388 lastlayer = l_layer[il];
389 increment++;
390 }
391 m_ioasym = (lout - lin) / increment;
392
393 // Get the vector of dE/dx values for all hits
394 if (m_isHitLevel) {
395 for (int ihit = 0; ihit < h_nhits; ++ihit) {
396
397 if (m_isDeadwire) continue;
398
399 h_lwire[ihit] = dedxTrack->getWireInLayer(ihit);
400 h_wire[ihit] = dedxTrack->getWire(ihit);
401 h_layer[ihit] = dedxTrack->getHitLayer(ihit);
402 h_path[ihit] = dedxTrack->getPath(ihit);
403 h_dedx[ihit] = dedxTrack->getDedx(ihit);
404 h_adcraw[ihit] = dedxTrack->getADCBaseCount(ihit);
405 h_adccorr[ihit] = dedxTrack->getADCCount(ihit);
406 h_doca[ihit] = dedxTrack->getDoca(ihit);
407 h_ndoca[ihit] = h_doca[ihit] / dedxTrack->getCellHalfWidth(ihit);
408 h_ndocaRS[ihit] = dedxTrack->getDocaRS(ihit) / dedxTrack->getCellHalfWidth(ihit);
409 h_enta[ihit] = dedxTrack->getEnta(ihit);
410 h_entaRS[ihit] = dedxTrack->getEntaRS(ihit);
411 h_driftT[ihit] = dedxTrack->getDriftT(ihit);
412 h_driftD[ihit] = dedxTrack->getDriftD(ihit);
413
414 // Get extra variable from tracking
415 if (m_isExtraVar) {
416 h_WeightPionHypo[ihit] = dedxTrack->getWeightPionHypo(ihit);
417 h_WeightKaonHypo[ihit] = dedxTrack->getWeightKaonHypo(ihit);
418 h_WeightProtonHypo[ihit] = dedxTrack->getWeightProtonHypo(ihit);
419 h_foundByTrackFinder[ihit] = dedxTrack->getFoundByTrackFinder(ihit);
420 }
421 // Get the calibration constants
422 h_cosCor[ihit] = m_DBCosineCor->getMean(h_layer[ihit], m_cosTheta);
423
424 h_facnladc[ihit] = dedxTrack->getNonLADCCorrection(ihit);
425 h_wireGain[ihit] = m_DBWireGains->getWireGain(h_wire[ihit]);
426 h_twodCor[ihit] = m_DB2DCell->getMean(h_layer[ihit], h_ndocaRS[ihit], h_entaRS[ihit]);
427 h_onedCor[ihit] = m_DB1DCell->getMean(h_layer[ihit], h_entaRS[ihit]);
428 }
429 }
430
431}
432
433//---------------------------------------------------------------------------------------
434void HitLevelInfoWriterModule::recalculateDedx(CDCDedxTrack* dedxTrack, map<int, vector<double>>& l_var,
435 double (&cdcChi)[Const::ChargedStable::c_SetSize])
436{
437 vector<double> newLayerHits;
438 double newLayerDe = 0, newLayerDx = 0;
439 int nhitscombined = 0; // number of hits combined per layer
440 int wirelongesthit = 0; // wire number of longest hit
441 double longesthit = 0; // path length of longest hit
442
443 for (int ihit = 0; ihit < h_nhits; ++ihit) {
444 int jadcbase = dedxTrack->getADCBaseCount(ihit);
445 int jLayer = dedxTrack->getHitLayer(ihit);
446 double jWire = dedxTrack->getWire(ihit);
447 double jNDocaRS = dedxTrack->getDocaRS(ihit) / dedxTrack->getCellHalfWidth(ihit);
448 double jEntaRS = dedxTrack->getEntaRS(ihit);
449 double jPath = dedxTrack->getPath(ihit);
450
451 double correction = dedxTrack->getScaleFactor() * dedxTrack->getRunGain() * dedxTrack->getTimeMean() *
452 dedxTrack->getCosineCorrection(ihit) * dedxTrack->getCosEdgeCorrection() *
453 dedxTrack->getTwoDCorrection(ihit) * dedxTrack->getOneDCorrection(ihit) *
454 dedxTrack->getNonLADCCorrection(ihit);
455 if (dedxTrack->getWireGain(ihit) > 0) correction *= dedxTrack->getWireGain(ihit); //also keep dead wire
456
457 if (m_isRelative) {
458 //get same base adc + rel correction factor
459 correction *= GetCorrection(jadcbase, jLayer, jWire, jNDocaRS, jEntaRS, m_cosTheta, m_injring, m_injtime);
460 if (!m_DBWireGains && dedxTrack->getWireGain(ihit) == 0) correction = 0;
461 } else {
462 //get modifed adc + abs correction factor
463 correction = GetCorrection(jadcbase, jLayer, jWire, jNDocaRS, jEntaRS, m_cosTheta, m_injring, m_injtime);
464 }
465
466 double newhitdedx = 1.0;
467
468 newhitdedx *= jadcbase * sqrt(1 - m_cosTheta * m_cosTheta) / jPath;
469
470 if (correction != 0) {
471 newhitdedx /= correction;
472 if (m_DBWireGains->getWireGain(jWire) != 0) {
473 newLayerDe += jadcbase / correction;
474 newLayerDx += jPath;
475 if (jPath > longesthit) {
476 longesthit = jPath;
477 wirelongesthit = jWire;
478 }
479 nhitscombined++;
480 }
481 } else newhitdedx = 0;
482
483 dedxTrack->setDedx(ihit, newhitdedx);
484
485 if (ihit + 1 < h_nhits && dedxTrack->getHitLayer(ihit + 1) == jLayer) {
486 continue;
487 } else {
488 if (newLayerDx != 0) {
489 double totalDistance = newLayerDx / sqrt(1 - m_cosTheta * m_cosTheta);
490 double newLayerDedx = newLayerDe / totalDistance ;
491 newLayerHits.push_back(newLayerDedx);
492 l_var[0].push_back(nhitscombined);
493 l_var[1].push_back(wirelongesthit);
494 l_var[2].push_back(jLayer);
495 l_var[3].push_back(totalDistance);
496 l_var[4].push_back(newLayerDedx);
497 }
498
499 newLayerDe = 0;
500 newLayerDx = 0;
501 nhitscombined = 0;
502 wirelongesthit = 0;
503 longesthit = 0;
504 }
505 }
506
507 // recalculate the truncated means
508 double dedxmean, dedxtrunc, dedxtruncNoSat, dedxerror;
509
510 calculateMeans(&dedxmean, &dedxtruncNoSat, &dedxerror, newLayerHits);
511
512 dedxtrunc = dedxtruncNoSat;
513
514 HadronCorrection(m_cosTheta, dedxtrunc);
515 m_mean = dedxmean;
516 m_trunc = dedxtrunc;
517 m_truncNoSat = dedxtruncNoSat;
518 m_error = dedxerror;
519
520 // save the PID information
521 saveChiValue(cdcChi, dedxTrack, m_trunc);
522}
523
524//---------------------------------------------------------------------------------------
525double HitLevelInfoWriterModule::GetCorrection(int& adc, int layer, int wireID, double doca, double enta,
526 double costheta, double ring, double time) const
527{
528 double correction = 1.0;
529 correction *= m_DBScaleFactor->getScaleFactor();
530 correction *= m_DBRunGain->getRunGain();
531 correction *= m_DB2DCell->getMean(layer, doca, enta);
532 correction *= m_DB1DCell->getMean(layer, enta);
533 correction *= m_DBCosineCor->getMean(layer, costheta);
534 correction *= m_DBInjectTime->getCorrection("mean", ring, time);
535 if (costheta <= -0.850 || costheta >= 0.950) correction *= m_DBCosEdgeCor->getMean(costheta);
536 if (m_DBWireGains->getWireGain(wireID) > 0) correction *= m_DBWireGains->getWireGain(wireID);
537
538 //last is only for abs constant
539 if (!m_isRelative) adc = m_DBNonlADC->getCorrectedADC(adc, layer);
540
541 return correction;
542}
543
544//---------------------------------------------------------------------------------------
545void HitLevelInfoWriterModule::HadronCorrection(double costheta, double& dedx) const
546{
547 dedx = D2I(costheta, I2D(costheta, 1.00) / 1.00 * dedx);
548}
549
550//---------------------------------------------------------------------------------------
551double HitLevelInfoWriterModule::D2I(const double cosTheta, const double D) const
552{
553 double absCosTheta = fabs(cosTheta);
554 double projection = pow(absCosTheta, m_hadronpars[3]) + m_hadronpars[2];
555 if (projection == 0) {
556 B2WARNING("Something wrong with dE/dx hadron constants!");
557 return D;
558 }
559
560 double chargeDensity = D / projection;
561 double numerator = 1 + m_hadronpars[0] * chargeDensity;
562 double denominator = 1 + m_hadronpars[1] * chargeDensity;
563
564 if (denominator == 0) {
565 B2WARNING("Something wrong with dE/dx hadron constants!");
566 return D;
567 }
568
569 double I = D * m_hadronpars[4] * numerator / denominator;
570 return I;
571}
572
573//---------------------------------------------------------------------------------------
574double HitLevelInfoWriterModule::I2D(const double cosTheta, const double I) const
575{
576 double absCosTheta = fabs(cosTheta);
577 double projection = pow(absCosTheta, m_hadronpars[3]) + m_hadronpars[2];
578
579 if (projection == 0 || m_hadronpars[4] == 0) {
580 B2WARNING("Something wrong with dE/dx hadron constants!");
581 return I;
582 }
583
584 double a = m_hadronpars[0] / projection;
585 double b = 1 - m_hadronpars[1] / projection * (I / m_hadronpars[4]);
586 double c = -1.0 * I / m_hadronpars[4];
587
588 if (b == 0 && a == 0) {
589 B2WARNING("both a and b coefficiants for hadron correction are 0");
590 return I;
591 }
592
593 double discr = b * b - 4.0 * a * c;
594 if (discr < 0) {
595 B2WARNING("negative discriminant; return uncorrectecd value");
596 return I;
597 }
598
599 double D = (a != 0) ? (-b + sqrt(discr)) / (2.0 * a) : -c / b;
600 if (D < 0) {
601 B2WARNING("D is less 0! will try another solution");
602 D = (a != 0) ? (-b - sqrt(discr)) / (2.0 * a) : -c / b;
603 if (D < 0) {
604 B2WARNING("D is still less 0! just return uncorrectecd value");
605 return I;
606 }
607 }
608
609 return D;
610}
611
612//---------------------------------------------------------------------------------------
613void HitLevelInfoWriterModule::calculateMeans(double* mean, double* truncMean, double* truncMeanErr,
614 const vector<double>& dedx) const
615{
616 // Calculate the truncated average by skipping the lowest & highest
617 // events in the array of dE/dx values
618 vector<double> sortedDedx = dedx;
619 sort(sortedDedx.begin(), sortedDedx.end());
620 sortedDedx.erase(remove(sortedDedx.begin(), sortedDedx.end(), 0), sortedDedx.end());
621 sortedDedx.shrink_to_fit();
622
623 double truncMeanTmp = 0.0, meanTmp = 0.0, sumOfSqs = 0.0;
624 int nValTrunc = 0;
625 const int numDedx = sortedDedx.size();
626
627 // add a factor of 0.51 here to make sure we are rounding appropriately...
628 const int lEdgeTrunc = int(numDedx * 0.05 + 0.51);
629 const int hEdgeTrunc = int(numDedx * (1 - 0.25) + 0.51);
630 for (int i = 0; i < numDedx; i++) {
631 meanTmp += sortedDedx[i];
632 if (i >= lEdgeTrunc and i < hEdgeTrunc) {
633 truncMeanTmp += sortedDedx[i];
634 sumOfSqs += sortedDedx[i] * sortedDedx[i];
635 nValTrunc++;
636 }
637 }
638
639 if (numDedx != 0) meanTmp /= numDedx;
640
641 if (nValTrunc != 0) truncMeanTmp /= nValTrunc;
642 else truncMeanTmp = meanTmp;
643
644 *mean = meanTmp;
645 *truncMean = truncMeanTmp;
646
647 if (nValTrunc > 1)
648 *truncMeanErr = sqrt(sumOfSqs / double(nValTrunc) - truncMeanTmp * truncMeanTmp) / double(nValTrunc - 1);
649 else *truncMeanErr = 0;
650
651}
652
653//---------------------------------------------------------------------------------------
655 double dedx) const
656{
657 // determine a chi value for each particle type
659 for (const Const::ChargedStable pdgIter : set) {
660
661 // determine the predicted mean and resolution
662 double mean = dedxTrack->getPmean(pdgIter.getIndex());
663 double sigma = dedxTrack->getPreso(pdgIter.getIndex());
664
665 // fill the chi value for this particle type
666 if (sigma != 0) chi[pdgIter.getIndex()] = ((dedx - mean) / (sigma));
667 }
668}
669
670//---------------------------------------------------------------------------------------
672{
673 for (int il = 0; il < 200; ++il) {
674 l_nhitscombined[il] = 0;
675 l_wirelongesthit[il] = 0;
676 l_layer[il] = 0;
677 l_path[il] = 0;
678 l_dedx[il] = 0;
679 l_cosCor[il] = 0;
680 }
681
682 if (m_isHitLevel) {
683 for (int ihit = 0; ihit < 200; ++ihit) {
684 h_lwire[ihit] = 0;
685 h_wire[ihit] = 0;
686 h_layer[ihit] = 0;
687 h_path[ihit] = 0;
688 h_dedx[ihit] = 0;
689 h_adcraw[ihit] = 0;
690 h_adccorr[ihit] = 0;
691 h_doca[ihit] = 0;
692 h_ndoca[ihit] = 0;
693 h_ndocaRS[ihit] = 0;
694 h_enta[ihit] = 0;
695 h_entaRS[ihit] = 0;
696 h_driftT[ihit] = 0;
697 // h_driftD[ihit] = 0;
698 h_facnladc[ihit] = 0;
699 h_cosCor[ihit] = 0;
700 h_wireGain[ihit] = 0;
701 h_twodCor[ihit] = 0;
702 h_onedCor[ihit] = 0;
703 h_WeightPionHypo[ihit] = 0;
704 h_WeightKaonHypo[ihit] = 0;
705 h_WeightProtonHypo[ihit] = 0;
706 h_foundByTrackFinder[ihit] = 0;
707 }
708 }
709}
710
711//---------------------------------------------------------------------------------------
713{
714 // register output root file
715 m_file.push_back(new TFile(filename.c_str(), "RECREATE"));
716 m_tree.push_back(new TTree("track", "dE/dx information"));
717
718 int i = m_tree.size() - 1;
719 m_tree[i]->SetDirectory(0);
720
721 // event level information (from emd)
722 m_tree[i]->Branch("exp", &m_expID, "exp/I");
723 m_tree[i]->Branch("run", &m_runID, "run/I");
724 m_tree[i]->Branch("event", &m_eventID, "event/I");
725
726 // track level information (from tfr)
727 m_tree[i]->Branch("dz", &m_dz, "dz/D");
728 m_tree[i]->Branch("dr", &m_dr, "dr/D");
729 m_tree[i]->Branch("chi2", &m_chi2, "chi2/D");
730 m_tree[i]->Branch("injtime", &m_injtime, "injtime/D");
731 m_tree[i]->Branch("isher", &m_injring, "isher/D");
732
733 if (m_isExtraVar) {
734 m_tree[i]->Branch("d0", &m_d0, "d0/D");
735 m_tree[i]->Branch("z0", &m_z0, "z0/D");
736 m_tree[i]->Branch("dphi", &m_dphi, "dphi/D");
737 m_tree[i]->Branch("vx0", &m_vx0, "vx0/D");
738 m_tree[i]->Branch("vy0", &m_vy0, "vy0/D");
739 m_tree[i]->Branch("vz0", &m_vz0, "vz0/D");
740 m_tree[i]->Branch("tanlambda", &m_tanlambda, "tanlambda/D");
741 m_tree[i]->Branch("phi0", &m_phi0, "phi0/D");
742 m_tree[i]->Branch("e1_9", &m_e1_9, "e1_9/D");
743 m_tree[i]->Branch("e9_21", &m_e9_21, "e9_21/D");
744 m_tree[i]->Branch("eclsnhits", &m_eclsnHits, "eclsnhits/D");
745 m_tree[i]->Branch("klmLayers", &m_klmLayers, "klmLayers/I");
746 }
747
748 // track level information (from cdt)
749 m_tree[i]->Branch("nCDChits", &m_nCDChits, "nCDChits/D");
750 m_tree[i]->Branch("inCDC", &m_inCDC, "inCDC/I");
751 m_tree[i]->Branch("track", &m_trackID, "track/I");
752 m_tree[i]->Branch("length", &m_length, "length/D");
753 m_tree[i]->Branch("charge", &m_charge, "charge/I");
754 m_tree[i]->Branch("costh", &m_cosTheta, "costh/D");
755 m_tree[i]->Branch("pF", &m_pCDC, "pF/D");
756 m_tree[i]->Branch("p", &m_p, "p/D");
757 m_tree[i]->Branch("pt", &m_pt, "pt/D");
758 m_tree[i]->Branch("ioasym", &m_ioasym, "ioasym/D");
759 m_tree[i]->Branch("phi", &m_phi, "phi/D");
760 m_tree[i]->Branch("pdg", &m_PDG, "pdg/D");
761 // track level dE/dx measurements
762 m_tree[i]->Branch("mean", &m_mean, "mean/D");
763 m_tree[i]->Branch("dedx", &m_trunc, "dedx/D");
764 m_tree[i]->Branch("dedxnosat", &m_truncNoSat, "dedxnosat/D");
765 m_tree[i]->Branch("dedxerr", &m_error, "dedxerr/D");
766
767 // other track level information
768 m_tree[i]->Branch("eop", &m_eop, "eop/D");
769 m_tree[i]->Branch("e", &m_e, "e/D");
770
771 // calibration constants
772 m_tree[i]->Branch("scale", &m_scale, "scale/D");
773 m_tree[i]->Branch("cosedgecor", &m_cosEdgeCor, "cosedgecor/D");
774 m_tree[i]->Branch("rungain", &m_runGain, "rungain/D");
775 m_tree[i]->Branch("timegain", &m_timeGain, "timegain/D");
776 m_tree[i]->Branch("timereso", &m_timeReso, "timereso/D");
777
778 // PID values
779 m_tree[i]->Branch("chiE", &m_chie, "chiE/D");
780 m_tree[i]->Branch("chiMu", &m_chimu, "chiMu/D");
781 m_tree[i]->Branch("chiPi", &m_chipi, "chiPi/D");
782 m_tree[i]->Branch("chiK", &m_chik, "chiK/D");
783 m_tree[i]->Branch("chiP", &m_chip, "chiP/D");
784 m_tree[i]->Branch("chiD", &m_chid, "chiD/D");
785
786 m_tree[i]->Branch("chiEOld", &m_chieOld, "chiEOld/D");
787 m_tree[i]->Branch("chiMuOld", &m_chimuOld, "chiMuOld/D");
788 m_tree[i]->Branch("chiPiOld", &m_chipiOld, "chiPiOld/D");
789 m_tree[i]->Branch("chiKOld", &m_chikOld, "chiKOld/D");
790 m_tree[i]->Branch("chiPOld", &m_chipOld, "chiPOld/D");
791 m_tree[i]->Branch("chiDOld", &m_chidOld, "chiDOld/D");
792
793 m_tree[i]->Branch("pmeanE", &m_pmeane, "pmeanE/D");
794 m_tree[i]->Branch("pmeanMu", &m_pmeanmu, "pmeanMu/D");
795 m_tree[i]->Branch("pmeanPi", &m_pmeanpi, "pmeanPi/D");
796 m_tree[i]->Branch("pmeanK", &m_pmeank, "pmeanK/D");
797 m_tree[i]->Branch("pmeanP", &m_pmeanp, "pmeanP/D");
798 m_tree[i]->Branch("pmeanD", &m_pmeand, "pmeanD/D");
799
800 m_tree[i]->Branch("presE", &m_prese, "presE/D");
801 m_tree[i]->Branch("presMu", &m_presmu, "presMu/D");
802 m_tree[i]->Branch("presPi", &m_prespi, "presPi/D");
803 m_tree[i]->Branch("presK", &m_presk, "presK/D");
804 m_tree[i]->Branch("presP", &m_presp, "presP/D");
805 m_tree[i]->Branch("presD", &m_presd, "presD/D");
806
807 // layer level information
808 m_tree[i]->Branch("lNHits", &l_nhits, "lNHits/I");
809 m_tree[i]->Branch("lNHitsUsed", &l_nhitsused, "lNHitsUsed/I");
810 m_tree[i]->Branch("lNHitsCombined", l_nhitscombined, "lNHitsCombined[lNHits]/I");
811 m_tree[i]->Branch("lWireLongestHit", l_wirelongesthit, "lWireLongestHit[lNHits]/I");
812 m_tree[i]->Branch("lLayer", l_layer, "lLayer[lNHits]/I");
813 m_tree[i]->Branch("lPath", l_path, "lPath[lNHits]/D");
814 m_tree[i]->Branch("lDedx", l_dedx, "lDedx[lNHits]/D");
815 m_tree[i]->Branch("lCosCor", l_cosCor, "lCosCor[lNHits]/D");
816
817 // hit level information
818 if (m_isHitLevel) {
819 m_tree[i]->Branch("hNHits", &h_nhits, "hNHits/I");
820 m_tree[i]->Branch("hLWire", h_lwire, "hLWire[hNHits]/I");
821 m_tree[i]->Branch("hWire", h_wire, "hWire[hNHits]/I");
822 m_tree[i]->Branch("hLayer", h_layer, "hLayer[hNHits]/I");
823 m_tree[i]->Branch("hPath", h_path, "hPath[hNHits]/D");
824 m_tree[i]->Branch("hDedx", h_dedx, "hDedx[hNHits]/D");
825 m_tree[i]->Branch("hADCRaw", h_adcraw, "hADCRaw[hNHits]/D");
826 m_tree[i]->Branch("hADCCorr", h_adccorr, "hADCCorr[hNHits]/D");
827 m_tree[i]->Branch("hDoca", h_doca, "hDoca[hNHits]/D");
828 m_tree[i]->Branch("hNDoca", h_ndoca, "hNDoca[hNHits]/D");
829 m_tree[i]->Branch("hNDocaRS", h_ndocaRS, "hNDocaRS[hNHits]/D");
830 m_tree[i]->Branch("hEnta", h_enta, "hEnta[hNHits]/D");
831 m_tree[i]->Branch("hEntaRS", h_entaRS, "hEntaRS[hNHits]/D");
832 m_tree[i]->Branch("hDriftT", h_driftT, "hDriftT[hNHits]/D");
833 m_tree[i]->Branch("hDriftD", h_driftD, "hDriftD[hNHits]/D");
834 m_tree[i]->Branch("hFacnlADC", h_facnladc, "hFacnlADC[hNHits]/D");
835 m_tree[i]->Branch("hCosCor", &h_cosCor, "hcosCor[hNHits]/D");
836 m_tree[i]->Branch("hWireGain", h_wireGain, "hWireGain[hNHits]/D");
837 m_tree[i]->Branch("hTwodcor", h_twodCor, "hTwodcor[hNHits]/D");
838 m_tree[i]->Branch("hOnedcor", h_onedCor, "hOnedcor[hNHits]/D");
839
840 if (m_isExtraVar) {
841 m_tree[i]->Branch("hWeightPionHypo", h_WeightPionHypo, "hWeightPionHypo[hNHits]/D");
842 m_tree[i]->Branch("hWeightKaonHypo", h_WeightKaonHypo, "hWeightKaonHypo[hNHits]/D");
843 m_tree[i]->Branch("hWeightProtonHypo", h_WeightProtonHypo, "hWeightProtonHypo[hNHits]/D");
844 m_tree[i]->Branch("hFoundByTrackFinder", h_foundByTrackFinder, "hFoundByTrackFinder[hNHits]/I");
845 }
846 }
847}
double R
typedef autogenerated by FFTW
Debug output for CDCDedxPID module.
double getTimeMean() const
Return the injection gain for this track.
int getHitLayer(int i) const
Return the (global) layer number for a hit.
int getADCBaseCount(int i) const
Return the base adcCount (no non-linearity) for this hit.
int getADCCount(int i) const
Return the adcCount for this hit.
double getDedxMean() const
Get the dE/dx mean for this track.
double getDriftD(int i) const
Return the drift distance for this hit.
int getLayer(int i) const
Return the (global) layer number for a layer hit.
double getDoca(int i) const
Return the distance of closest approach to the sense wire for this hit.
double getDedx() const
Get dE/dx truncated mean for this track.
double getPath(int i) const
Return the path length through the cell for this hit.
int getWireInLayer(int i) const
Return the sensor ID for this hit: wire number in the layer.
double getPmean(int i) const
Return the PID (predicted mean) value.
double getCosineCorrection(int i) const
Return the cosine correction for this hit.
double getOneDCorrection(int i) const
Return the 1D correction for this hit.
int getWireLongestHit(int i) const
Return the wire number of the longest hit per layer.
double getWeightKaonHypo(int i) const
Return the max weights from KalmanFitterInfo using kaon hypothesis.
void setDedx(double mean)
Set the dE/dx truncated average for this track.
double getEntaRS(int i) const
Return rescaled enta value for cell height=width assumption.
double getLayerDedx(int i) const
Return the total dE/dx for this layer.
double getCellHalfWidth(int i) const
Return the half-width of the CDC cell.
double getPDG() const
Get the identity of the particle.
double getTwoDCorrection(int i) const
Return the 2D correction for this hit.
double getPreso(int i) const
Return the PID (predicted reso) value.
int getNLayerHits() const
Return the number of layer hits for this track.
int getFoundByTrackFinder(int i) const
Return the TrackFinder which added the given hit to track.
int getDriftT(int i) const
Return the drift time for this hit.
double getCosTheta() const
Return cos(theta) for this track.
double getInjectionRing() const
Return cos(theta) for this track.
int getCharge() const
Return the charge for this track.
double getWireGain(int i) const
Return the wire gain for this hit.
double getWeightPionHypo(int i) const
Return the max weights from KalmanFitterInfo using pion hypothesis.
int getWire(int i) const
Return the sensor ID for this hit: wire number for CDC (0-14336)
double getDocaRS(int i) const
Return rescaled doca value for cell height=width assumption.
double getCosEdgeCorrection() const
Return the cosine correction for this track.
double getRunGain() const
Return the run gain for this track.
double getNonLADCCorrection(int i) const
Return the factor introduce for adcCount (non-linearity) correction.
double getDedxError() const
Get the error on the dE/dx truncated mean for this track.
double getWeightProtonHypo(int i) const
Return the max weights from KalmanFitterInfo using proton hypothesis.
double getEnta(int i) const
Return the entrance angle in the CDC cell for this hit.
double getInjectionTime() const
Return cos(theta) for this track.
double trackID() const
Return the track ID.
double getLayerPath(int i) const
Return the distance travelled in this layer.
double getChi(int i) const
Return the PID (chi) value.
double getScaleFactor() const
Return the scale factor for this track.
double getLength() const
Return the total path length for this track.
int getNHitsCombined(int i) const
Return the number of hits combined per layer.
double getDedxNoSat() const
Get dE/dx truncated mean without the saturation correction for this track.
int size() const
Return the number of hits for this track.
double getNLayerHitsUsed() const
Return the number of hits used to determine the truncated mean.
double getMomentum() const
Return the track momentum valid in the CDC.
Provides a type-safe way to pass members of the chargedStableSet set.
Definition Const.h:589
static const unsigned int c_SetSize
Number of elements (for use in array bounds etc.)
Definition Const.h:615
A set of ParticleType objects, with defined order.
Definition Const.h:517
static const ChargedStable muon
muon particle
Definition Const.h:660
static const ParticleSet chargedStableSet
set of charged stable particles
Definition Const.h:618
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
Class for accessing objects in the database.
Definition DBObjPtr.h:21
ECL cluster data.
Definition ECLCluster.h:27
bool hasHypothesis(EHypothesisBit bitmask) const
Return if specific hypothesis bit is set.
Definition ECLCluster.h:360
double getE1oE9() const
Return E1/E9 (shower shape variable).
Definition ECLCluster.h:286
double getE9oE21() const
Return E9/E21 (shower shape variable).
Definition ECLCluster.h:289
double getEnergy(EHypothesisBit hypothesis) const
Return Energy (GeV).
Definition ECLCluster.cc:23
double getNumberOfCrystals() const
Return number of a crystals in a shower (sum of weights).
Definition ECLCluster.h:304
@ c_nPhotons
CR is split into n photons (N1)
Definition ECLCluster.h:41
double h_ndoca[kMaxHits]
normalized distance of closest approach
double h_driftD[kMaxHits]
drift distance
double h_onedCor[kMaxHits]
calibration 1D cleanup correction
int l_nhits
the total number of layer hits for this Track
double m_cosEdgeCor
calibration cosine edge correction
double m_presk
pred reso for kaon hypothesis
double m_phi0
Angle of the transverse momentum in the r-phi plane.
double m_d0
Signed distance to the POCA in the r-phi plane.
double m_presp
pred reso for proton hypothesis
double D2I(const double cosTheta, const double D) const
hadron saturation parameterization part 2
int l_wirelongesthit[kMaxHits]
the wire number of longest hit in this layer
double m_ioasym
asymmetry in increasing vs decreasing layer numbers per track
std::vector< double > m_hadronpars
hadron saturation parameters
double m_p
momentum from tracking
double m_e9_21
ratio of energies of the central 3x3 crystal vs 5x5 crystals
int m_trackID
ID number of the Track.
DBObjPtr< CDCDedxRunGain > m_DBRunGain
Run gain DB object.
double m_vx0
X coordinate of track POCA to origin.
virtual ~HitLevelInfoWriterModule()
Destructor.
StoreArray< KLMCluster > m_klmClusters
Required array of input KLMClusters.
double m_pmeand
pred mean for deuteron hypothesis
virtual void initialize() override
Initialize the module.
DBObjPtr< CDCDedxHadronCor > m_DBHadronCor
hadron saturation parameters
double h_enta[kMaxHits]
entrance angle
double m_chipOld
chi value for proton hypothesis
std::string m_strOutputBaseName
Base name for the output ROOT files.
virtual void event() override
This method is called for each event.
double h_WeightKaonHypo[kMaxHits]
weight for kaon hypothesis from KalmanFitterInfo
DBObjPtr< CDCDedxCosineEdge > m_DBCosEdgeCor
cosine edge calibration
HitLevelInfoWriterModule()
Default constructor.
double m_pmeanp
pred mean for proton hypothesis
double I2D(const double cosTheta, const double I) const
hadron saturation parameterization part 1
double h_WeightProtonHypo[kMaxHits]
weight for proton hypothesis from KalmanFitterInfo
double m_chid
modified chi value for deuteron hypothesis
double m_chikOld
chi value for kaon hypothesis
double m_dz
vertex or POCA in case of tracks z in respect to IPs
bool m_isHitLevel
Flag to switch on/off hit level info.
DBObjPtr< CDCDedxADCNonLinearity > m_DBNonlADC
hadron saturation non linearity
double m_pmeank
pred mean for kaon hypothesis
double h_facnladc[kMaxHits]
calibration hit gain
DBObjPtr< CDCDedx1DCell > m_DB1DCell
1D correction DB object
void saveChiValue(double(&chi)[Const::ChargedStable::c_SetSize], CDCDedxTrack *dedxTrack, double dedx) const
for all particles, save chi values into 'chi' chi array of chi values to be modified
StoreArray< TrackFitResult > m_trackFitResults
Required array of input TrackFitResults.
int h_nhits
the number of good hits for this Track
virtual void terminate() override
End of the event processing.
double h_adccorr[kMaxHits]
charge per hit corr by nonlinear ADC
std::vector< TTree * > m_tree
output ROOT trees
double m_presmu
pred reso for muon hypothesis
double m_tanlambda
Slope of the track in the r-z plane.
int h_layer[kMaxHits]
layer number
double h_adcraw[kMaxHits]
charge per hit
int m_expID
experiment in which this Track was found
void fillTrack(const TrackFitResult *fitResult)
Fill the TTree with the information from the track fit.
double h_path[kMaxHits]
path length in cell
double m_chimu
modified chi value for muon hypothesis
double m_vy0
Y coordinate of track POCA to origin.
double m_chik
modified chi value for kaon hypothesis
void recalculateDedx(CDCDedxTrack *dedxTrack, std::map< int, std::vector< double > > &l_var, double(&cdcChi)[Const::ChargedStable::c_SetSize])
Function to recalculate the dedx with latest constants.
double m_chi2
chi^2 from track fit
void HadronCorrection(double costheta, double &dedx) const
Function to apply the hadron correction.
void calculateMeans(double *mean, double *truncatedMean, double *truncatedMeanErr, const std::vector< double > &dedx) const
Save arithmetic and truncated mean for the 'dedx' values.
double h_wireGain[kMaxHits]
calibration hit gain
double m_theta
cos(theta) for the track
double h_twodCor[kMaxHits]
calibration 2D correction
double m_chie
modified chi value for electron hypothesis
DBObjPtr< CDCDedx2DCell > m_DB2DCell
2D correction DB object
double m_runGain
calibration run gain
int m_eventID
event in which this Track was found
double m_timeGain
calibration injection time gain
void bookOutput(std::string filename)
Create the output TFiles and TTrees.
double m_chidOld
chi value for deuteron hypothesis
bool m_isExtraVar
Flag to switch on/off extra level info and some available w/ release/5 only.
void clearEntries()
Clear the arrays before filling an event.
std::vector< TFile * > m_file
output ROOT files
double m_dr
track d0 relative to IP
double m_prese
pred reso for electron hypothesis
DBObjPtr< CDCDedxCosineCor > m_DBCosineCor
Electron saturation correction DB object.
int h_foundByTrackFinder[kMaxHits]
the 'found by track finder' flag for the given hit
double m_timeReso
calibration injection time reso
StoreArray< Track > m_tracks
Required array of input Tracks.
double l_dedx[kMaxHits]
dE/dx for this layer
double l_path[kMaxHits]
distance travelled in this layer
int l_layer[kMaxHits]
layer number
double m_chipiOld
chi value for pion hypothesis
double m_chip
modified chi value for proton hypothesis
bool m_isCorrection
Flag to switch on/off corrections.
double m_pmeane
pred mean for electron hypothesis
double m_cosTheta
cos(theta) for the track
int l_nhitsused
the total number of layer hits used for this Track
StoreArray< ECLCluster > m_eclClusters
Required array of input ECLClusters.
double h_dedx[kMaxHits]
charge per path length
double m_chieOld
chi value for electron hypothesis
double h_driftT[kMaxHits]
drift time
double m_dphi
POCA in degrees in respect to IP.
DBObjPtr< CDCDedxWireGain > m_DBWireGains
Wire gain DB object.
double m_eop
energy over momentum in the calorimeter
bool m_isDeadwire
write only active wires
double m_presd
pred reso for deuteron hypothesis
double m_scale
calibration scale factor
double m_klmLayers
number of klm layers with hits
int l_nhitscombined[kMaxHits]
the number of hits combined this layer
double m_error
standard deviation of the truncated mean
double m_pCDC
momentum valid in CDC
double m_truncNoSat
dE/dx averaged, truncated mean, with corrections (not hadron)
double m_length
total path length of the Track
double m_injring
HER injection status.
double m_nCDChits
Number of CDC hits associated to the track.
double h_ndocaRS[kMaxHits]
normalized +RS distance of closest approach
double m_vz0
Z coordinate of track POCA to origin.
double h_cosCor[kMaxHits]
calibration cosine correction
double m_trunc
dE/dx averaged, truncated mean, with corrections
double GetCorrection(int &adc, int layer, int wireID, double doca, double enta, double costheta, double ring, double time) const
Function to get the correction factor.
int m_inCDC
frack is CDC acceptance or not
double h_doca[kMaxHits]
distance of closest approach
double l_cosCor[kMaxHits]
calibration cosine correction
double m_chimuOld
chi value for muon hypothesis
double m_z0
z coordinate of the POCA
int m_charge
the charge for this Track
double h_WeightPionHypo[kMaxHits]
weight for pion hypothesis from KalmanFitterInfo
int m_runID
run in which this Track was found
double m_e
energy in the calorimeter
StoreArray< CDCDedxTrack > m_dedxTracks
Required array of CDCDedxTracks.
double m_injtime
time since last injection in micro seconds
double m_pmeanmu
pred mean for muon hypothesis
double h_entaRS[kMaxHits]
normalized + RS distance of entrance angle
bool m_isRelative
Flag to switch on/off relative constants.
double m_prespi
pred reso for pion hypothesis
DBObjPtr< CDCDedxInjectionTime > m_DBInjectTime
time gain/reso DB object
int h_wire[kMaxHits]
sense wire ID
DBObjPtr< CDCDedxScaleFactor > m_DBScaleFactor
Scale factor to make electrons ~1.
int h_lwire[kMaxHits]
sense wire within layer
double m_chipi
modified chi value for pion hypothesis
std::vector< std::string > m_strParticleList
Vector of ParticleLists to write out.
double m_pt
transverse momentum from tracking
double m_pmeanpi
pred mean for pion hypothesis
void fillDedx(CDCDedxTrack *dedxTrack)
Fill the TTree with the information from a CDCDedxTrack object.
double m_e1_9
ratio of energies of the central 1 crystal vs 3x3 crystals
unsigned short getNHits() const
Get the total Number of CDC hits in the fit.
KLM cluster data.
Definition KLMCluster.h:29
int getLayers() const
Get number of layers with hits.
Definition KLMCluster.h:67
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
Module()
Constructor.
Definition Module.cc:30
Class to collect log likelihoods from TOP, ARICH, dEdx, ECL and KLM aimed for output to mdst includes...
Class to store reconstructed particles.
Definition Particle.h:76
static const ReferenceFrame & GetCurrent()
Get current rest frame.
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
T * getRelated(const std::string &name="", const std::string &namedRelation="") const
Get the object to or from which this object has a relation.
FROM * getRelatedFrom(const std::string &name="", const std::string &namedRelation="") const
Get the object from which this object has a relation.
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.
double getD0() const
Getter for d0.
double getTanLambda() const
Getter for tanLambda.
double getZ0() const
Getter for z0.
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;.
UncertainHelix getUncertainHelix() const
Conversion to framework Uncertain Helix (i.e., with covariance).
double getPhi0() const
Getter for phi0.
Class that bundles various TrackFitResults.
Definition Track.h:25
This class represents an ideal helix in perigee parameterization including the covariance matrix of t...
double passiveMoveBy(const ROOT::Math::XYZVector &by)
Moves origin of the coordinate system (passive transformation) by the given vector.
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:559
#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.