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