CbmRoot
Loading...
Searching...
No Matches
CbmTofFindTracks.cxx
Go to the documentation of this file.
1/* Copyright (C) 2015-2021 PI-UHd, GSI
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Norbert Herrmann [committer], Pierre-Alain Loizeau */
4
5// -------------------------------------------------------------------------
6// ----- CbmTofFindTracks source file -----
7// ----- Created 25/04/15 by N. Herrmann -----
8// ----- initially following CbmTrdFindTracks -----
9// -------------------------------------------------------------------------
10
11#include "CbmTofFindTracks.h"
12
13#include "CbmEvent.h"
14#include "CbmMatch.h"
15#include "CbmTofAddress.h" // in cbmdata/tof
16#include "CbmTofCalibrator.h"
17#include "CbmTofCell.h" // in tof/TofData
18#include "CbmTofCreateDigiPar.h" // in tof/TofTools
19#include "CbmTofDetectorId_v12b.h" // in cbmdata/tof
20#include "CbmTofDetectorId_v14a.h" // in cbmdata/tof
21#include "CbmTofDetectorId_v21a.h" // in cbmdata/tof
22#include "CbmTofDigiBdfPar.h" // in tof/TofParam
23#include "CbmTofDigiPar.h" // in tof/TofParam
24#include "CbmTofGeoHandler.h" // in tof/TofTools
25#include "CbmTofHit.h"
26#include "CbmTofTrackFinder.h"
27#include "CbmTofTrackFinderNN.h"
28#include "CbmTofTracklet.h"
29#include "CbmTofTrackletTools.h"
30#include "CbmVertex.h"
31#include "FairRootManager.h"
32#include "FairRunAna.h"
33#include "FairRuntimeDb.h"
34#include "TClonesArray.h"
35#include "TDirectory.h"
36#include "TF1.h"
37#include "TFile.h"
38#include "TFitResult.h"
39#include "TGeoManager.h"
40#include "TH1.h"
41#include "TH2.h"
42#include "TH3.h"
43#include "TMath.h"
44#include "TProfile.h"
45#include "TROOT.h"
46#include "TRandom.h"
47#include "TSpectrum.h"
48#include "TString.h"
49
50#include <Logger.h>
51
52#include <iostream>
53#include <vector>
54
55using std::vector;
56
57//const Int_t DetMask = 0x3FFFFF; // check for consistency with v14a geometry
58const Int_t DetMask = 0x1FFFFF; // check for consistency with v21a geometry
59static int iTS = 0;
60
62
64
65// ----- Default constructor -------------------------------------------
67{
68 if (!fInstance) fInstance = this;
69}
70// -------------------------------------------------------------------------
71
72
73// ----- Standard constructor ------------------------------------------
74CbmTofFindTracks::CbmTofFindTracks(const char* name, const char* /*title*/, CbmTofTrackFinder* finder)
75 : FairTask(name)
76 , fFinder(finder)
77 , fTrackletTools(NULL)
78 , fTofCalibrator(NULL)
79 , fEventsColl(NULL)
80 , fTofHitArrayIn(NULL)
81 , fTofMatchArrayIn(NULL)
82 , fTofHitArray(NULL)
83 , fTofHitIndexArray()
84 , fTofHitArrayOut(NULL)
85 , fTofUHitArrayOut(NULL)
86 , fTrackArray(NULL)
87 , fTrackArrayOut(nullptr)
88 , fTofUHitArray(NULL)
89 , fMinNofHits(-1)
90 , fNofTracks(-1)
91 , fNTofStations(-1)
92 , fNReqStations(-1)
93 , fInspectEvent(kTRUE)
94 , fStationType()
95 , fStationHMul()
96 , fRpcAddr()
97 , fMapStationRpcId()
98 , fMapRpcIdParInd()
99 , fvToff()
100 , fvXoff()
101 , fvYoff()
102 , fvZoff()
103 , fvTsig()
104 , fvXsig()
105 , fvYsig()
106 , fvZsig()
107 , fhTrklMul(NULL)
108 , fhTrklChi2(NULL)
109 , fhAllHitsStation(NULL)
110 , fhAllHitsSmTypes(NULL)
111 , fhUsedHitsStation(NULL)
112 , fhTrackingTimeNhits(NULL)
113 , fhTrklMulNhits(NULL)
114 , fhTrklMulMaxMM(NULL)
115 , fhTrklMul3D(NULL)
116 , fhTrklHMul(NULL)
117 , fhTrklZ0xHMul(NULL)
118 , fhTrklZ0yHMul(NULL)
119 , fhTrklTxHMul(NULL)
120 , fhTrklTyHMul(NULL)
121 , fhTrklTyTx(NULL)
122 , fhTrklTtHMul(NULL)
123 , fhTrklVelHMul(NULL)
124 , fhTrklT0HMul(NULL)
125 , fhTrklT0Mul(NULL)
126 , fhTrklDT0SmMis(NULL)
127 , fhTrklDT0StMis2(NULL)
128 , fhTrklXY0_0(NULL)
129 , fhTrklXY0_1(NULL)
130 , fhTrklXY0_2(NULL)
131 , fhTrklX0_TX(NULL)
132 , fhTrklY0_TX(NULL)
133 , fhTrklX0_TY(NULL)
134 , fhTrklY0_TY(NULL)
135 , vhPullX()
136 , vhPullY()
137 , vhPullZ()
138 , vhPullT()
139 , vhPullTB()
140 , vhTrefRms()
141 , vhFitDT0()
142 , vhFitT0Err()
143 , vhFitTt()
144 , vhFitTtErr()
145 , vhFitDTMean()
146 , vhFitDTMeanErr()
147 , vhResidualTBWalk()
148 , vhResidualYWalk()
149 , vhXY_AllTracks()
150 , vhXY_AllStations()
151 , vhXY_AllFitStations()
152 , vhXY_MissedStation()
153 , vhXY_DX()
154 , vhXY_DY()
155 , vhXY_DT()
156 , vhXY_TOT()
157 , vhXY_CSZ()
158 , vhUDXDY_DT()
159 , vhUCDXDY_DT()
160 , fhVTXNorm(NULL)
161 , fhVTX_XY0(NULL)
162 , fhVTX_DT0_Norm(NULL)
163 , fiStationStatus(0)
164 , fOutHstFileName("")
165 , fCalParFileName("")
166 , fCalOutFileName("./tofFindTracks.hst.root")
167 , fCalParFile(NULL)
168 , fhPullT_Smt(NULL)
169 , fhPullT_Smt_Off(NULL)
170 , fhPullX_Smt(NULL)
171 , fhPullX_Smt_Off(NULL)
172 , fhPullY_Smt(NULL)
173 , fhPullY_Smt_Off(NULL)
174 , fhPullZ_Smt(NULL)
175 , fhPullZ_Smt_Off(NULL)
176 , fhPullT_Smt_Width(NULL)
177 , fhPullX_Smt_Width(NULL)
178 , fhPullY_Smt_Width(NULL)
179 , fhPullZ_Smt_Width(NULL)
180 , fhTOff_Smt(NULL)
181 , fhTOff_Smt_Off(NULL)
182 , fhDeltaTt_Smt(NULL)
183 , fhDeltaTc_Smt(NULL)
184 , fhTOff_HMul2(NULL)
185 , fiCorMode(0)
186 , fiBeamCounter(-1)
187 , fiStationMaxHMul(1000)
188 , fTtTarg(30.)
189 , fTtLight(0.)
190 , fdTOffScal(1.)
191 , fVTXNorm(0.)
192 , fVTX_T(0.)
193 , fVTX_X(0.)
194 , fVTX_Y(0.)
195 , fVTX_Z(0.)
196 , fT0MAX(0.5)
197 , fiEvent(0)
198 , fGeoHandler(new CbmTofGeoHandler())
199 , fTofId(NULL)
200 , fDigiPar(NULL)
201 , fDigiBdfPar(NULL)
202 , fSIGT(0.1)
203 , fSIGX(1.)
204 , fSIGY(1.)
205 , fSIGZ(1.)
206 , fbUseSigCalib(kTRUE)
207 , fdRefVelMean(0.)
208 , fdRefDVel(1.E7)
209 , fdR0Lim(0.)
210 , fdTtMin(0.)
211 , fStart()
212 , fStop()
213 , fdTrackingTime(0.)
214 , fdBeamMomentumLab(0.)
215 , fbRemoveSignalPropagationTime(kFALSE)
216 , fiBeamMaxHMul(1000)
217 , fiCalOpt((int) 0)
218{
219 if (!fInstance) fInstance = this;
220}
221// -------------------------------------------------------------------------
222
223
224// ----- Destructor ----------------------------------------------------
226{
227 if (fInstance == this) fInstance = 0;
228 fTrackArray->Delete();
229}
230// -------------------------------------------------------------------------
231
232
233// ----- Public method Init (abstract in base class) --------------------
235{
236
237 // Check for Track finder
238 if (!fFinder) {
239 LOG(warning) << "-W- CbmTofFindTracks::Init: No track finder selected!";
240 return kERROR;
241 }
242
243 fTrackletTools = new CbmTofTrackletTools(); // initialize tools
244
245 // Get and check FairRootManager
246 FairRootManager* ioman = FairRootManager::Instance();
247 if (!ioman) {
248 LOG(error) << "-E- CbmTofFindTracks::Init: RootManager not instantiated!";
249 return kFATAL;
250 }
251
252 ioman->InitSink();
253
254 fEventsColl = dynamic_cast<TClonesArray*>(ioman->GetObject("Event"));
255 if (!fEventsColl) fEventsColl = dynamic_cast<TClonesArray*>(ioman->GetObject("CbmEvent"));
256 if (!fEventsColl) {
257 LOG(info) << "CbmEvent not found in input file, assume eventwise input";
258 }
259
260 fTofHitArray = new TClonesArray("CbmTofHit");
261
262
263 // Get TOF hit Array
264 fTofHitArrayIn = (TClonesArray*) ioman->GetObject("TofHit");
265 if (!fTofHitArrayIn) {
266 LOG(fatal) << "-W- CbmTofFindTracks::Init: No TofHit array!";
267 return kERROR;
268 }
269
270 // Get TOF DigiMatch Array
271 fTofMatchArrayIn = (TClonesArray*) ioman->GetObject("TofHitCalDigiMatch");
272 if (!fTofMatchArrayIn) {
273 LOG(fatal) << "CbmTofFindTracks::Init: No TofDigiMatch array!";
274 return kERROR;
275 }
276
277 // Create and register output TofTrack array
278 fTrackArray = new TClonesArray("CbmTofTracklet", 100);
279 fTofUHitArray = new TClonesArray("CbmTofHit", 100);
280 //fTrackArray->BypassStreamer(kTRUE); //needed?
281 //ioman->Register("TofTracklets", "TOF", fTrackArray, kFALSE); //FIXME
282 if (fEventsColl) {
283 fTrackArrayOut = new TClonesArray("CbmTofTracklet", 100);
284 fTofHitArrayOut = new TClonesArray("CbmTofHit", 100);
285 fTofUHitArrayOut = new TClonesArray("CbmTofHit", 100);
286 ioman->Register("TofTracklets", "TOF", fTrackArrayOut, IsOutputBranchPersistent("TofTracklets"));
287 ioman->Register("TofCalHit", "TOF", fTofHitArrayOut, IsOutputBranchPersistent("TofCalHit"));
288 LOG(info) << "-W- TofTracklets array registered in EventsColl mode";
289 ioman->Register("TofUHit", "TOF", fTofUHitArrayOut, IsOutputBranchPersistent("TofTracklets"));
290 }
291 else {
292 ioman->Register("TofTracklets", "TOF", fTrackArray, IsOutputBranchPersistent("TofTracklets"));
293 LOG(info) << "-W- CbmTofFindTracks::Init:TofTracklets array registered";
294 LOG(info) << "Register TofCalHit at " << fTofHitArray;
295 ioman->Register("TofCalHit", "TOF", fTofHitArray, IsOutputBranchPersistent("TofCalHit"));
296
297 // Create and register TofUHit array for unused Hits
298 ioman->Register("TofUHit", "TOF", fTofUHitArray, kFALSE);
299 }
300 // Call the Init method of the track finder
301 fFinder->Init();
302
303 if (fOutHstFileName == "") {
304 fOutHstFileName = "./FindTofTracks.hst.root";
305 }
306 LOG(info) << "CbmTofFindTracks::Init: Hst Output filename = " << fOutHstFileName;
307
308 if (kFALSE == InitParameters()) return kFATAL;
309
310 // default parameters
311 // if (fMinNofHits < 1) fMinNofHits=1;
312
313 //fill RpcId - map
314 Bool_t bBeamCounter = kFALSE;
315 Int_t iRpc = 0;
316 fRpcAddr.resize(0);
317 for (Int_t iCell = 0; iCell < fDigiPar->GetNrOfModules(); iCell++) {
318 Int_t iCellId = fDigiPar->GetCellId(iCell);
319 Int_t iCh = fTofId->GetCell(iCellId);
320 if (0 == iCh) {
321 LOG(info) << Form("Init found at Ind %d, %lu Rpc with Addr 0x%08x, TSR %d%d%d ", iRpc, fRpcAddr.size(), iCellId,
322 fTofId->GetSMType(iCellId), fTofId->GetSModule(iCellId), fTofId->GetCounter(iCellId));
323 if (fTofId->GetSMType(iCellId) == 5) {
324 bBeamCounter = kTRUE;
325 LOG(info) << "Found beam counter in setup! at RpcInd " << iRpc << ", Addr.size " << fRpcAddr.size();
326 }
327 fMapRpcIdParInd[iCellId] = iRpc;
328 fRpcAddr.push_back(iCellId);
329 iRpc++;
330 }
331 }
332
333 LOG(debug) << "Initialize fStationHMul to size " << fNTofStations + 1;
334 fStationHMul.resize(fNTofStations + 1);
335 LOG(debug) << "Initialize fStationHMul to size " << fNTofStations + 1;
336
338
340
341 if (fiCalOpt > 0) {
343 if (fTofCalibrator->Init() != kSUCCESS) return kFATAL;
344 if (bBeamCounter) {
345 if (fiBeamCounter > -1) fTofCalibrator->SetBeam(bBeamCounter);
347 LOG(info) << "Set CbmTofCalibrator::R0Lim to " << fdR0Lim;
348 }
349 }
350
351 LOG(info) << Form("BeamCounter to be used in tracking: 0x%08x", fiBeamCounter);
352
353 return kSUCCESS;
354}
355// -------------------------------------------------------------------------
356/************************************************************************************/
358{
359 UInt_t NSt = fMapRpcIdParInd.size();
360 fvToff.resize(NSt);
361 fvXoff.resize(NSt);
362 fvYoff.resize(NSt);
363 fvZoff.resize(NSt);
364 fvTsig.resize(NSt);
365 fvXsig.resize(NSt);
366 fvYsig.resize(NSt);
367 fvZsig.resize(NSt);
368 for (uint i = 0; i < NSt; i++) {
369 fvToff[i] = 0.;
370 fvXoff[i] = 0.;
371 fvYoff[i] = 0.;
372 fvZoff[i] = 0.;
373 fvTsig[i] = fSIGT;
374 fvXsig[i] = fSIGX;
375 fvYsig[i] = fSIGY;
376 fvZsig[i] = fSIGZ;
377 }
378
379 if (fCalParFileName.IsNull()) return kTRUE;
380
382 TFile* oldFile = gFile;
383 TDirectory* oldDir = gDirectory;
384
385 fCalParFile = new TFile(fCalParFileName, "");
386 if (NULL == fCalParFile) {
387 LOG(error) << "CbmTofFindTracks::LoadCalParameter: "
388 << "file " << fCalParFileName << " does not exist!";
389 return kTRUE;
390 }
391
392 LOG(info) << "CbmTofFindTracks::LoadCalParameter: "
393 << " read from file " << fCalParFileName;
394
395 TH1D* fhtmp = (TH1D*) gDirectory->FindObjectAny(Form("hPullT_Smt_Off"));
396 TH1D* fhtmpX = (TH1D*) gDirectory->FindObjectAny(Form("hPullX_Smt_Off"));
397 TH1D* fhtmpY = (TH1D*) gDirectory->FindObjectAny(Form("hPullY_Smt_Off"));
398 TH1D* fhtmpZ = (TH1D*) gDirectory->FindObjectAny(Form("hPullZ_Smt_Off"));
399 TH1D* fhtmpWT = (TH1D*) gDirectory->FindObjectAny(Form("hPullT_Smt_Width"));
400 TH1D* fhtmpWX = (TH1D*) gDirectory->FindObjectAny(Form("hPullX_Smt_Width"));
401 TH1D* fhtmpWY = (TH1D*) gDirectory->FindObjectAny(Form("hPullY_Smt_Width"));
402 TH1D* fhtmpWZ = (TH1D*) gDirectory->FindObjectAny(Form("hPullZ_Smt_Width"));
403
404 gROOT->cd();
405 if (NULL == fhtmp) {
406 LOG(info) << Form("CbmTofFindTracks::LoadCalParameter: hPullT_Smt_Off") << " not found. ";
407 }
408 else {
409 fhPullT_Smt_Off = (TH1D*) fhtmp->Clone();
410 for (UInt_t iSt = 0; iSt < NSt; iSt++)
411 fvToff[iSt] = fhPullT_Smt_Off->GetBinContent(iSt + 1);
412 }
413
414 if (NULL == fhtmpX) {
415 LOG(info) << Form("CbmTofFindTracks::LoadCalParameter: hPullX_Smt_Off") << " not found. ";
416 }
417 else {
418 fhPullX_Smt_Off = (TH1D*) fhtmpX->Clone();
419 for (UInt_t iSt = 0; iSt < NSt; iSt++)
420 fvXoff[iSt] = fhPullX_Smt_Off->GetBinContent(iSt + 1);
421 }
422
423 if (NULL == fhtmpY) {
424 LOG(info) << Form("CbmTofFindTracks::LoadCalParameter: hPullY_Smt_Off") << " not found. ";
425 }
426 else {
427 fhPullY_Smt_Off = (TH1D*) fhtmpY->Clone();
428 for (UInt_t iSt = 0; iSt < NSt; iSt++)
429 fvYoff[iSt] = fhPullY_Smt_Off->GetBinContent(iSt + 1);
430 }
431
432 if (NULL == fhtmpZ) {
433 LOG(info) << Form("CbmTofFindTracks::LoadCalParameter: hPullZ_Smt_Off") << " not found. ";
434 }
435 else {
436 fhPullZ_Smt_Off = (TH1D*) fhtmpZ->Clone();
437 for (UInt_t iSt = 0; iSt < NSt; iSt++)
438 fvZoff[iSt] = fhPullZ_Smt_Off->GetBinContent(iSt + 1);
439 }
440
441 if (NULL == fhtmpWT) {
442 LOG(info) << Form("CbmTofFindTracks::LoadCalParameter: hPullT_Smt_Width") << " not found. ";
443 }
444 else {
445 if (fbUseSigCalib) {
446 fhPullT_Smt_Width = (TH1D*) fhtmpWT->Clone();
447 for (UInt_t iSt = 0; iSt < NSt; iSt++) {
448 fvTsig[iSt] = fhPullT_Smt_Width->GetBinContent(iSt + 1);
449 if (fvTsig[iSt] == 0) {
450 LOG(warning) << "Invalid Tsig for station " << iSt;
451 fvTsig[iSt] = fSIGT;
452 }
453 }
454 }
455 }
456
457 if (NULL == fhtmpWX) {
458 LOG(info) << Form("CbmTofFindTracks::LoadCalParameter: hPullX_Smt_Width") << " not found. ";
459 }
460 else {
461 if (fbUseSigCalib) {
462 fhPullX_Smt_Width = (TH1D*) fhtmpWX->Clone();
463 for (UInt_t iSt = 0; iSt < NSt; iSt++) {
464 fvXsig[iSt] = fhPullX_Smt_Width->GetBinContent(iSt + 1);
465 if (fvXsig[iSt] == 0) {
466 LOG(warning) << "Invalid Xsig for station " << iSt;
467 fvXsig[iSt] = fSIGX;
468 }
469 }
470 }
471 }
472
473 if (NULL == fhtmpWY) {
474 LOG(info) << Form("CbmTofFindTracks::LoadCalParameter: hPullY_Smt_Width") << " not found. ";
475 }
476 else {
477 if (fbUseSigCalib) {
478 fhPullY_Smt_Width = (TH1D*) fhtmpWY->Clone();
479 for (UInt_t iSt = 0; iSt < NSt; iSt++) {
480 fvYsig[iSt] = fhPullY_Smt_Width->GetBinContent(iSt + 1);
481 if (fvYsig[iSt] == 0) {
482 LOG(warning) << "Invalid Ysig for station " << iSt;
483 fvYsig[iSt] = fSIGY;
484 }
485 }
486 }
487 }
488
489 if (NULL == fhtmpWZ) {
490 LOG(info) << Form("CbmTofFindTracks::LoadCalParameter: hPullZ_Smt_Width") << " not found. ";
491 }
492 else {
493 if (fbUseSigCalib) {
494 fhPullZ_Smt_Width = (TH1D*) fhtmpWZ->Clone();
495 for (UInt_t iSt = 0; iSt < NSt; iSt++) {
496 fvZsig[iSt] = fhPullZ_Smt_Width->GetBinContent(iSt + 1);
497 if (fvZsig[iSt] == 0) {
498 LOG(warning) << "Invalid Zsig for station " << iSt;
499 fvZsig[iSt] = fSIGZ;
500 }
501 }
502 }
503 }
504
505 fCalParFile->Close();
506
507 Double_t nSmt = fMapRpcIdParInd.size();
508
509 if (NULL == fhPullT_Smt_Off) { // provide default TOffset histogram
511 new TH1F(Form("hPullT_Smt_Off"), Form("Tracklet PullT vs RpcInd ; RpcInd ; #DeltaT (ns)"), nSmt, 0, nSmt);
512 }
513 // Initialize Parameter
514 if (fiCorMode <= -3) { // hidden option, FIXME, disabled
515 for (Int_t iDet = 0; iDet < nSmt; iDet++) {
516 std::map<Int_t, Int_t>::iterator it;
517 //it = fMapRpcIdParInd.find(iDet)
518 Int_t iMap = 0;
519 for (it = fMapRpcIdParInd.begin(); it != fMapRpcIdParInd.end(); it++) {
520 iMap++;
521 if (it->second == iDet) break;
522 }
523 LOG(debug1) << Form(" iDet %d -> iUniqueId ? 0x%08x, 0x%08x ", iDet, it->first, it->second);
524 Int_t iUniqueId = it->first;
525 CbmTofCell* fChannelInfo = fDigiPar->GetCell(iUniqueId);
526 if (NULL != fChannelInfo) {
527 Double_t dVal = 0.; // FIXME numeric constant in code, default for cosmic
528 dVal = fhPullT_Smt_Off->GetBinContent(iDet + 1, dVal);
529 if (dVal == 0) {
530 if (fiBeamCounter != iUniqueId) dVal = fChannelInfo->GetZ() * fTtTarg; // use calibration target value
531 fhPullT_Smt_Off->SetBinContent(iDet + 1, dVal);
532 LOG(info) << Form("PrimInit det 0x%08x at %d, GloInd %d, z=%f to Tt %6.4f, %6.4f with TOff %6.2f", iUniqueId,
533 iDet, iMap, fChannelInfo->GetZ(), fTtTarg, fdTOffScal, dVal);
534 }
535 else {
536 if (fdTOffScal != 0.) fhPullT_Smt_Off->SetBinContent(iDet + 1, dVal * fdTOffScal);
537 }
538 LOG(info) << Form("ReInit det 0x%08x at %d, GloInd %d, z=%f to Tt %6.4f, %6.4f with TOff %6.2f", iUniqueId,
539 iDet, iMap, fChannelInfo->GetZ(), fTtTarg, fdTOffScal, dVal);
540 }
541 }
542 }
543 if (NULL == fhPullT_Smt_Width) { // provide default TWidth histogram
545 new TH1F(Form("hPullT_Smt_Width"), Form("Tracklet ResiT Width vs RpcInd ; RpcInd ; RMS(T) (ns)"), nSmt, 0, nSmt);
546
547 // Initialize Parameter
548 for (Int_t iDet = 0; iDet < nSmt; iDet++) {
549 fhPullT_Smt_Width->SetBinContent(iDet + 1, fSIGT);
550 }
551 }
552
553 LOG(info) << "CbmTofFindTracks::LoadCalParameter: fhPullT_Smt_Off at " << fhPullT_Smt_Off;
554
555 if (NULL == fhPullX_Smt_Off) // provide default XOffset histogram
557 new TH1F(Form("hPullX_Smt_Off"), Form("Tracklet ResiX vs RpcInd ; RpcInd ; #DeltaX (cm)"), nSmt, 0, nSmt);
558 if (NULL == fhPullX_Smt_Width) {
560 new TH1F(Form("hPullX_Smt_Width"), Form("Tracklet ResiX Width vs RpcInd ; RpcInd ; RMS(X) (cm)"), nSmt, 0, nSmt);
561 // Initialize Parameter
562 for (Int_t iDet = 0; iDet < nSmt; iDet++) {
563 fhPullX_Smt_Width->SetBinContent(iDet + 1, fSIGX);
564 }
565 }
566
567 if (NULL == fhPullY_Smt_Off) // provide default YOffset histogram
569 new TH1F(Form("hPullY_Smt_Off"), Form("Tracklet ResiY vs RpcInd ; RpcInd ; #DeltaY (cm)"), nSmt, 0, nSmt);
570 if (NULL == fhPullY_Smt_Width) {
572 new TH1F(Form("hPullY_Smt_Width"), Form("Tracklet ResiY Width vs RpcInd ; RpcInd ; RMS(Y) (cm)"), nSmt, 0, nSmt);
573 // Initialize Parameter
574 for (Int_t iDet = 0; iDet < nSmt; iDet++) {
575 fhPullY_Smt_Width->SetBinContent(iDet + 1, fSIGY);
576 }
577 }
578
579 if (NULL == fhPullZ_Smt_Off) // provide default ZOffset histogram
581 new TH1F(Form("hPullZ_Smt_Off"), Form("Tracklet ResiZ vs RpcInd ; RpcInd ; #DeltaZ (cm)"), nSmt, 0, nSmt);
582 if (NULL == fhPullZ_Smt_Width) {
584 new TH1F(Form("hPullZ_Smt_Width"), Form("Tracklet ResiZ Width vs RpcInd ; RpcInd ; RMS(Z) (cm)"), nSmt, 0, nSmt);
585 // Initialize Parameter
586 for (Int_t iDet = 0; iDet < nSmt; iDet++) {
587 fhPullZ_Smt_Width->SetBinContent(iDet + 1, fSIGZ);
588 }
589 }
590
592 gFile = oldFile;
593 gDirectory = oldDir;
594
595 return kTRUE;
596}
597//-------------------------------------------------------------------------------------------------
599{
600 // Initialize the TOF GeoHandler
601 Bool_t isSimulation = kFALSE;
602 Int_t iGeoVersion = fGeoHandler->Init(isSimulation);
603 if (k12b > iGeoVersion) {
604 LOG(error) << "CbmTofFindTracks::InitParameters => Only compatible with "
605 "geometries after v12b !!!";
606 return kFALSE;
607 }
608
609 LOG(info) << "CbmTofFindTracks::InitParameters: GeoVersion " << iGeoVersion;
610
611 switch (iGeoVersion) {
612 case k12b: fTofId = new CbmTofDetectorId_v12b(); break;
613 case k14a: fTofId = new CbmTofDetectorId_v14a(); break;
614 case k21a: fTofId = new CbmTofDetectorId_v21a(); break;
615 default: LOG(fatal) << "CbmTofFindTracks::InitParameters: Invalid Detector ID " << iGeoVersion;
616 }
617
618 // create digitization parameters from geometry file
619 CbmTofCreateDigiPar* tofDigiPar = new CbmTofCreateDigiPar("TOF Digi Producer", "TOF task");
620 LOG(info) << "Create DigiPar ";
621 tofDigiPar->Init();
622
623 return kTRUE;
624}
625// ----- SetParContainers -------------------------------------------------
627{
628 FairRunAna* ana = FairRunAna::Instance();
629 FairRuntimeDb* rtdb = ana->GetRuntimeDb();
630 // rtdb->getContainer("CbmGeoPassivePar");
631 // rtdb->getContainer("CbmGeoStsPar");
632 // rtdb->getContainer("CbmGeoTofPar");
633 rtdb->getContainer("FairBaseParSet");
634 // rtdb->getContainer("CbmGeoPassivePar");
635 // rtdb->getContainer("CbmGeoStsPar");
636 // rtdb->getContainer("CbmGeoRichPar");
637 rtdb->getContainer("CbmGeoTofPar");
638 // rtdb->getContainer("CbmFieldPar");
639 fDigiPar = (CbmTofDigiPar*) (rtdb->getContainer("CbmTofDigiPar"));
640
641 LOG(info) << " CbmTofFindTracks::SetParContainers found " << fDigiPar->GetNrOfModules() << " cells ";
642
643 fDigiBdfPar = (CbmTofDigiBdfPar*) (rtdb->getContainer("CbmTofDigiBdfPar"));
644}
645// -------------------------------------------------------------------------
646
648{
649 if (fiCorMode < 0) return kTRUE;
650
651 LOG(info) << Form("CbmTofFindTracks::WriteHistos: %s, mode = %d", fCalOutFileName.Data(), fiCorMode);
652
653 // Write histogramms to the file
655 TFile* oldFile = gFile;
656 TDirectory* oldDir = gDirectory;
657 TFile* fHist = new TFile(fCalOutFileName, "RECREATE");
658 fHist->cd();
659 const Double_t RMSmin = 0.03; // in ns
660
661 switch (fiCorMode) {
662 case 0: {
663 TProfile* htmp = fhPullT_Smt->ProfileX();
664 TH1D* htmp1D = htmp->ProjectionX();
665 TProfile* hTOff = fhTOff_HMul2->ProfileX();
666 TH1D* hTOff1D = hTOff->ProjectionX();
667
668 Double_t nx = htmp1D->GetNbinsX();
669 for (Int_t ix = 1; ix < nx; ix++) {
670 Double_t dVal = 0;
671 if (fhPullT_Smt_Off != NULL) {
672 dVal = fhPullT_Smt_Off->GetBinContent(ix + 1);
673 }
674 else {
675 fhPullT_Smt_Off = htmp1D;
676 }
677 TH1D* hTOff1DY = fhTOff_HMul2->ProjectionY(Form("_py%d", ix), ix + 1, ix + 1, "");
678 Double_t dFMean = 0.;
679 if (hTOff1DY->GetEntries() > 100) {
680 //Double_t dMean=hTOff1DY->GetMean();
681 Int_t iBmax = hTOff1DY->GetMaximumBin();
682 TAxis* xaxis = hTOff1DY->GetXaxis();
683 Double_t dMean = xaxis->GetBinCenter(iBmax); //X-value of bin with maximal content
684 Double_t dLim = 1000.; //1.5*hTOff1DY->GetRMS();
685 TFitResultPtr fRes = hTOff1DY->Fit("gaus", "S", "", dMean - dLim, dMean + dLim);
686 dFMean = fRes->Parameter(1);
687 }
688 dVal -= dFMean;
689 LOG(info) << "Init TOff " << ix << ": Old " << fhPullT_Smt_Off->GetBinContent(ix + 1) << ", Cnts "
690 << hTOff1D->GetBinContent(ix + 1) << ", FitMean " << dFMean << " -> " << dVal;
691 fhPullT_Smt_Off->SetBinContent(ix + 1, dVal);
692 }
693 }
694
695 break;
696
697 case 1: // correct mean deviation from fit (Pull)
698 {
699 TProfile* htmp = fhPullT_Smt->ProfileX();
700 TH1D* htmp1D = htmp->ProjectionX();
701 TProfile* hTOff = fhTOff_Smt->ProfileX();
702 TH1D* hTOff1D = hTOff->ProjectionX();
703
704 if (fhPullT_Smt_Off != NULL) {
705 Double_t nx = htmp1D->GetNbinsX();
706 for (Int_t ix = 0; ix < nx; ix++) {
707 Double_t dVal = fhPullT_Smt_Off->GetBinContent(ix + 1);
708 dVal -= htmp1D->GetBinContent(ix + 1);
709
710 LOG(debug1) << "Update hPullT_Smt_Off " << ix << ": " << fhPullT_Smt_Off->GetBinContent(ix + 1) << " + "
711 << htmp1D->GetBinContent(ix + 1) << " + " << hTOff1D->GetBinContent(ix + 1) << " -> " << dVal;
712 fhPullT_Smt_Off->SetBinContent(ix + 1, dVal);
713 }
714 }
715 else {
716 LOG(warning) << "CbmTofFindTracks::WriteHistos: fhPullT_Smt_Off not found ";
717 }
718 }
719
720 break;
721
722 case 2: // correct deviation from DeltaTt=0 expectation
723 {
724 TProfile* htmp = fhPullT_Smt->ProfileX();
725 TH1D* htmp1D = htmp->ProjectionX();
726 TProfile* hTOff = fhTOff_Smt->ProfileX();
727 TH1D* hTOff1D = hTOff->ProjectionX();
728
729 if (fhPullT_Smt_Off != NULL) {
730 Double_t nx = htmp1D->GetNbinsX();
731 for (Int_t ix = 0; ix < nx; ix++) {
732 Double_t dVal = fhPullT_Smt_Off->GetBinContent(ix + 1);
733 dVal -= hTOff1D->GetBinContent(ix + 1);
734 TH1D* hTOff1DY = fhTOff_Smt->ProjectionY(Form("_py%d", ix), ix + 1, ix + 1, "");
735 Double_t dFMean = 0.;
736 if (hTOff1DY->GetEntries() > 100) {
737 //Double_t dMean=hTOff1DY->GetMean();
738 Int_t iBmax = hTOff1DY->GetMaximumBin();
739 TAxis* xaxis = hTOff1DY->GetXaxis();
740 Double_t dMean = xaxis->GetBinCenter(iBmax); //X-value of bin with maximal content
741 Double_t dLim = 1.5 * hTOff1DY->GetRMS();
742 TFitResultPtr fRes = hTOff1DY->Fit("gaus", "SQM", "", dMean - dLim, dMean + dLim);
743 Int_t iFitStatus = fRes;
744 //if (iFitStatus == 0) { // check validity of fit
745 dFMean = fRes->Parameter(1);
746 dVal += hTOff1D->GetBinContent(ix + 1); //revert default correction
747 dVal -= dFMean;
748 //}
749 LOG(info) << "Update hPullT_Smt_Off Ind " << ix << ", stat " << iFitStatus << ": Old "
750 << fhPullT_Smt_Off->GetBinContent(ix + 1) << ", Pull " << htmp1D->GetBinContent(ix + 1)
751 << ", Dev@Peak " << hTOff1D->GetBinContent(ix + 1) << ", FitMean " << dFMean << " -> " << dVal;
752 }
753 else {
754 LOG(debug1) << "Update hPullT_Smt_Off " << ix << ": insufficient counts: " << hTOff1DY->GetEntries();
755 }
756 fhPullT_Smt_Off->SetBinContent(ix + 1, dVal);
757 }
758 }
759 else {
760 LOG(warning) << "CbmTofFindTracks::WriteHistos: fhPullT_Smt_Off not found ";
761 }
762 } break;
763
764 case 3: // correct Time Offset from PullT, extract width
765 {
766 TProfile* htmp = fhPullT_Smt->ProfileX();
767 TH1D* htmp1D = htmp->ProjectionX();
768 if (fhPullT_Smt_Off != NULL) {
769 Double_t nx = htmp1D->GetNbinsX();
770 for (Int_t ix = 0; ix < nx; ix++) {
771 Double_t dVal = fhPullT_Smt_Off->GetBinContent(ix + 1); //Current value
772 Double_t dCor = htmp->GetBinContent(ix + 1);
773 Double_t dRMS = htmp->GetBinError(ix + 1);
774 TH1D* hpy = fhPullT_Smt->ProjectionY(Form("%s_py%d", fhPullT_Smt->GetName(), ix), ix + 1, ix + 1);
775 if (hpy->GetEntries() > 50.) {
776 Int_t iBmax = hpy->GetMaximumBin();
777 TAxis* xaxis = hpy->GetXaxis();
778 Double_t dMean = xaxis->GetBinCenter(iBmax); //X-value of bin with maximal content
779 dRMS = TMath::Abs(hpy->GetRMS());
780 Double_t dLim = 1.5 * dRMS;
781 Double_t dNorm = hpy->GetBinContent(iBmax);
782 LOG(info) << "Fit3 " << hpy->GetName()
783 << Form(", %f with %f, %f, %f ", hpy->GetEntries(), dNorm, dMean, dLim);
784 if (dNorm > 10) {
785 TFitResultPtr fRes = hpy->Fit("gaus", "SQM", "", dMean - dLim, dMean + dLim);
786 //TF1* mgaus=new TF1("mgaus","gaus", dMean - dLim, dMean + dLim);
787 //mgaus->SetParameters(dNorm,dMean,dLim*0.5);
788 //TFitResultPtr fRes = hpy->Fit("mgaus", "SQM", "", dMean - dLim, dMean + dLim);
789 // see https://root-forum.cern.ch/t/tfitresultptr-not-valid-check/35944/4
790 if (gMinuit->fCstatu.Contains("OK") || gMinuit->fCstatu.Contains("CONVERGED")) {
791 TF1* fg = hpy->GetFunction("gaus");
792 if (fg == NULL) {
793 LOG(fatal) << "No associated gaus function for " << hpy->GetName();
794 continue;
795 }
796 //Double_t dFMean = fRes->Parameter(1);
797 Double_t dFMean = fg->GetParameter(1);
798 dCor = dFMean; // update offset
799 Double_t dFMeanError = fg->GetParError(1);
800 LOG(info) << "Update hPullT_Smt_Off3 Ind " << ix << Form(", 0x%08x: ", fRpcAddr[ix])
801 << fhPullT_Smt_Off->GetBinContent(ix + 1) << " + " << dFMean << ", Err " << dFMeanError
802 << " -> " << dVal - dCor << ", Width " << dRMS << ", Chi2 " << fg->GetChisquare();
803 if (dFMeanError < 0.05) { // FIXME: hardwired constant
804 if (dRMS < RMSmin) dRMS = RMSmin;
805 if (dRMS > fSIGT * 3.0) dRMS = fSIGT * 3.;
806 }
807 }
808 else {
809 LOG(info) << " Fit of " << hpy->GetName() << " failed with " << gMinuit->fCstatu;
810 }
811 }
812 else {
813 LOG(info) << "Fit3: Too few entries for fit ofhisto " << hpy->GetName() << ": " << dNorm;
814 }
815 }
816 else {
817 LOG(info) << "Update hPullT_Smt_Off " << ix << ": insufficient counts: " << hpy->GetEntries();
818 }
819
820 if (fRpcAddr[ix] != fiBeamCounter) // don't correct beam counter time
821 fhPullT_Smt_Off->SetBinContent(ix + 1, dVal - dCor);
822 else
823 LOG(info) << "No Off3 correction for beam counter at index " << ix;
824
825 fhPullT_Smt_Width->SetBinContent(ix + 1, dRMS);
826 } //ix loop end
827 }
828 else {
829 LOG(warning) << "CbmTofFindTracks::WriteHistos: fhPullT_Smt_Off not found ";
830 }
831 }
832
833 break;
834
835 case 4: // correct mean deviation from fit (Pull), extract width for x direction
836 {
837 TProfile* htmp = fhPullX_Smt->ProfileX();
838 TH1D* htmp1D = htmp->ProjectionX();
839
840 if (fhPullX_Smt_Off != NULL) {
841 Double_t nx = htmp1D->GetNbinsX();
842 for (Int_t ix = 0; ix < nx; ix++) {
843 int iSmType = CbmTofAddress::GetSmType(fRpcAddr[ix] & DetMask);
844 if (iSmType == 8) continue; // skip pad counters
845 TH1D* hpy = fhPullX_Smt->ProjectionY("_py", ix + 1, ix + 1);
846 Double_t dVal = fhPullX_Smt_Off->GetBinContent(ix + 1);
847 //dVal -= htmp1D->GetBinContent(ix + 1);
848 if (hpy->GetEntries() > 100.) {
849 // Fit gaussian
850 Double_t dFMean = hpy->GetBinCenter(hpy->GetMaximumBin());
851 Double_t dFLim = 0.5; // CAUTION, fixed numeric value
852 Double_t dBinSize = hpy->GetBinWidth(1);
853 dFLim = TMath::Max(dFLim, 5. * dBinSize);
854 TFitResultPtr fRes = hpy->Fit("gaus", "SQM", "", dFMean - dFLim, dFMean + dFLim);
855 dVal -= fRes->Parameter(1);
856 Double_t dRMS = fRes->Parameter(2);
857 LOG(debug) << "PeakFit at " << dFMean << ", lim " << dFLim << " : mean " << fRes->Parameter(1) << ", width "
858 << dRMS;
859
860 dRMS = TMath::Abs(hpy->GetRMS());
861 if (dRMS < fSIGX * 0.5) dRMS = fSIGX * 0.5;
862 if (dRMS > fSIGX * 3.0) dRMS = fSIGX * 3.;
863
864 // limit maximal shift in X, for larger values, change geometry file
865 const double dMaxShift = 4.;
866 if (dVal < -dMaxShift) dVal = -dMaxShift;
867 if (dVal > dMaxShift) dVal = dMaxShift;
868 LOG(info) << "Update hPullX_Smt_Off " << ix << ": " << fhPullX_Smt_Off->GetBinContent(ix + 1) << " + "
869 << htmp1D->GetBinContent(ix + 1) << ", " << fRes->Parameter(1) << " -> " << dVal << ", Width "
870 << dRMS;
871 if (fRpcAddr[ix] != fiBeamCounter) // don't correct beam counter position
872 fhPullX_Smt_Off->SetBinContent(ix + 1, dVal);
873 fhPullX_Smt_Width->SetBinContent(ix + 1, dRMS);
874 }
875 }
876 }
877 else {
878 LOG(warning) << "CbmTofFindTracks::WriteHistos: fhPullX_Smt_Off not found ";
879 }
880 }
881
882 break;
883
884 case 5: // correct mean deviation from fit (Pull), extract width for Y direction
885 {
886 TProfile* htmp = fhPullY_Smt->ProfileX();
887 TH1D* htmp1D = htmp->ProjectionX();
888
889 if (fhPullY_Smt_Off != NULL) {
890 Double_t nx = htmp1D->GetNbinsX();
891 for (Int_t ix = 0; ix < nx; ix++) {
892 int iSmType = CbmTofAddress::GetSmType(fRpcAddr[ix] & DetMask);
893 if (iSmType == 8) continue; // skip pad counters
894 Double_t dVal = fhPullY_Smt_Off->GetBinContent(ix + 1);
895 //dVal -= htmp1D->GetBinContent(ix + 1);
896 // Fit gaussian
897 TH1D* hpy = fhPullY_Smt->ProjectionY("_py", ix + 1, ix + 1);
898 if (hpy->GetEntries() > 100.) {
899 Double_t dFMean = hpy->GetBinCenter(hpy->GetMaximumBin());
900 Double_t dFLim = 2.; // CAUTION, fixed numeric value
901 Double_t dBinSize = hpy->GetBinWidth(1);
902 dFLim = TMath::Max(dFLim, 5. * dBinSize);
903 TFitResultPtr fRes = hpy->Fit("gaus", "SQM", "", dFMean - dFLim, dFMean + dFLim);
904 dVal -= fRes->Parameter(1);
905 Double_t dRMS = fRes->Parameter(2);
906 LOG(debug) << "PeakFit at " << dFMean << ", lim " << dFLim << " : mean " << fRes->Parameter(1) << ", width "
907 << dRMS;
908 if (fRpcAddr[ix] != fiBeamCounter) // don't correct beam counter position
909 fhPullY_Smt_Off->SetBinContent(ix + 1, dVal);
910
911 dRMS = TMath::Abs(hpy->GetRMS());
912 if (dRMS < fSIGY * 0.5) dRMS = 0.5 * fSIGY;
913 if (dRMS > fSIGY * 3.0) dRMS = fSIGY * 3.;
914 fhPullY_Smt_Width->SetBinContent(ix + 1, dRMS);
915
916 LOG(debug1) << "Update hPullY_Smt_Off " << ix << ": " << fhPullY_Smt_Off->GetBinContent(ix + 1) << " + "
917 << htmp1D->GetBinContent(ix + 1) << " -> " << dVal << ", Width " << dRMS;
918 }
919 }
920 }
921 else {
922 LOG(warning) << "CbmTofFindTracks::WriteHistos: fhPullY_Smt_Off not found ";
923 }
924
925 }
926
927 break;
928
929 case 6: // correct mean deviation from fit (Pull), extract width
930 {
931 TProfile* htmp = fhPullZ_Smt->ProfileX();
932 TH1D* htmp1D = htmp->ProjectionX();
933
934 if (fhPullZ_Smt_Off != NULL) {
935 Double_t nx = htmp1D->GetNbinsX();
936 for (Int_t ix = 0; ix < nx; ix++) {
937 Double_t dVal = fhPullZ_Smt_Off->GetBinContent(ix + 1);
938 dVal -= htmp1D->GetBinContent(ix + 1);
939 fhPullZ_Smt_Off->SetBinContent(ix + 1, dVal);
940
941 TH1D* hpy = fhPullZ_Smt->ProjectionY("_py", ix + 1, ix + 1);
942 if (hpy->GetEntries() > 100.) {
943 Double_t dRMS = TMath::Abs(hpy->GetRMS());
944
945 LOG(debug1) << "Update hPullZ_Smt_Off " << ix << ": " << fhPullZ_Smt_Off->GetBinContent(ix + 1) << " + "
946 << htmp1D->GetBinContent(ix + 1) << " -> " << dVal << ", Width " << dRMS;
947 if (dRMS < 1.5) dRMS = 1.5;
948 fhPullZ_Smt_Width->SetBinContent(ix + 1, dRMS);
949 }
950 }
951 }
952 else {
953 LOG(warning) << "CbmTofFindTracks::WriteHistos: fhPullZ_Smt_Off not found ";
954 }
955
956 } break;
957
958 case 7: // extract residual widthes in T, X, Y, Z
959 {
960 for (Int_t iStation = 0; iStation < static_cast<Int_t>(fMapRpcIdParInd.size()); iStation++) {
961 TH1D* hResidualT = fhPullT_Smt->ProjectionY("_py", iStation + 1, iStation + 1);
962 TH1D* hResidualX = fhPullX_Smt->ProjectionY("_py", iStation + 1, iStation + 1);
963 TH1D* hResidualY = fhPullY_Smt->ProjectionY("_py", iStation + 1, iStation + 1);
964 TH1D* hResidualZ = fhPullZ_Smt->ProjectionY("_py", iStation + 1, iStation + 1);
965
966 if (hResidualT->GetEntries() > 100.) {
967 Double_t dRMS = TMath::Abs(hResidualT->GetRMS());
968
969 if (dRMS < RMSmin) dRMS = RMSmin;
970 if (dRMS > 3. * fSIGT) dRMS = 3. * fSIGT;
971
972 fhPullT_Smt_Width->SetBinContent(iStation + 1, dRMS);
973 }
974
975 if (hResidualX->GetEntries() > 100.) {
976 Double_t dRMS = TMath::Abs(hResidualX->GetRMS());
977
978 if (dRMS < 0.5 * fSIGX) dRMS = 0.5 * fSIGX;
979 if (dRMS > 3. * fSIGX) dRMS = 3. * fSIGX;
980
981 fhPullX_Smt_Width->SetBinContent(iStation + 1, dRMS);
982 }
983
984 if (hResidualY->GetEntries() > 100.) {
985 Double_t dRMS = TMath::Abs(hResidualY->GetRMS());
986
987 if (dRMS < 0.5 * fSIGY) dRMS = 0.5 * fSIGY;
988 if (dRMS > 3. * fSIGY) dRMS = 3. * fSIGY;
989
990 fhPullY_Smt_Width->SetBinContent(iStation + 1, dRMS);
991 }
992
993 if (hResidualZ->GetEntries() > 100.) {
994 Double_t dRMS = TMath::Abs(hResidualZ->GetRMS());
995
996 if (dRMS < 1.5) dRMS = 1.5;
997
998 fhPullZ_Smt_Width->SetBinContent(iStation + 1, dRMS);
999 }
1000 }
1001 } break;
1002
1003 case 10: //correct mean deviation from TB - histo of station 0
1004 case 11:
1005 case 12:
1006 case 13:
1007 case 14:
1008 case 15:
1009 case 16:
1010 case 17:
1011 case 18:
1012 case 19: {
1013 Int_t iSt = fiCorMode % 10;
1014 TString hname = Form("hPull%s_Station_%d", "TB", iSt);
1015 TH1* h1 = (TH1*) gROOT->FindObjectAny(hname);
1016 if (h1->GetEntries() > 100) {
1017 Double_t dFMean = h1->GetMean();
1018 Double_t dFLim = 2.5 * h1->GetRMS();
1019 TFitResultPtr fRes = h1->Fit("gaus", "SQM", "", dFMean - dFLim, dFMean + dFLim);
1020 Double_t dDOff = fRes->Parameter(1);
1021 Double_t dSig = fRes->Parameter(2);
1022 Int_t iRpcInd = fMapRpcIdParInd[fMapStationRpcId[iSt]];
1023 Double_t dVal = fhPullT_Smt_Off->GetBinContent(iRpcInd + 1);
1024 dVal -= dDOff;
1025 LOG(info) << "Update hPullT_Smt_OffP Ind " << iSt << ", Ind " << iRpcInd << ": "
1026 << fhPullT_Smt_Off->GetBinContent(iRpcInd + 1) << " - " << dDOff << " -> " << dVal << ", Width "
1027 << dSig;
1028 fhPullT_Smt_Off->SetBinContent(iRpcInd + 1, dVal);
1029 if (dSig < fSIGT * 0.5) dSig = 0.5 * fSIGT;
1030 if (dSig > fSIGT * 3.0) dSig = fSIGT * 3.;
1031 fhPullT_Smt_Width->SetBinContent(iRpcInd + 1, dSig);
1032 }
1033 else {
1034 LOG(info) << "CbmTofFindTracks::WriteHistos: Too few entries in histo " << hname;
1035 }
1036 } break;
1037
1038 case 82:
1039 case 81:
1040 case 80: {
1041 Int_t iSel = fiCorMode % 10;
1042 LOG(info) << "Update time offsets with Detector Doublets " << iSel;
1043 int iO[3] = {0, 5, 10}; // 0 - 1 - 2 - layers
1044 switch (iSel) {
1045 case 1: iO[2] = 25; break; // 0 - 1 - big modules
1046 case 2:
1047 iO[0] = 5;
1048 iO[1] = 25;
1049 iO[2] = 12;
1050 break; // 1 - big - 2
1051 default:;
1052 }
1053 const size_t N = 3;
1054 double dTshift[N] = {3 * 0.};
1055 for (int iLoc = 0; iLoc < 5; iLoc++) { // loop over rpcs in module
1056 if (EvalDoublets(iLoc, iLoc + iO[1], iLoc + iO[2], dTshift)) { // returns vector of shifts
1057 double dTMeanShift = (dTshift[0] + dTshift[1] + dTshift[2]) / 3.;
1058 for (int i = 0; i < 3; i++) { // apply time shifts
1059 int iStation = GetStationOfAddr(fDigiBdfPar->GetDetUId(iLoc + iO[i]));
1060 if (fiStationStatus[iStation] > 0) continue; // do not modify
1061 int ix = fMapRpcIdParInd[fDigiBdfPar->GetDetUId(iLoc + iO[i])]; // convert BDF to Geo counting
1062 LOG(info) << "UpdateDT0 bdf ch " << iLoc + iO[i] << ", geo ch " << ix << " by " << Form("%f", dTshift[i]);
1063 fhPullT_Smt_Off->SetBinContent(ix + 1, fhPullT_Smt_Off->GetBinContent(ix + 1) + dTshift[i] - dTMeanShift);
1064 }
1065 }
1066 else {
1067 iO[2]++; // try neighbor
1068 if (EvalDoublets(iLoc, iLoc + iO[1], iLoc + iO[2], dTshift)) { // returns vector of shifts
1069 double dTMeanShift = (dTshift[0] + dTshift[1] + dTshift[2]) / 3.;
1070 for (int i = 0; i < 3; i++) {
1071 int iStation = GetStationOfAddr(fDigiBdfPar->GetDetUId(iLoc + iO[i]));
1072 if (fiStationStatus[iStation] > 0) continue; // do not modify
1073 int ix = fMapRpcIdParInd[fDigiBdfPar->GetDetUId(iLoc + iO[i])]; // convert BDF to Geo counting
1074 LOG(info) << "UpdateDT+ bdf ch " << iLoc + iO[i] << ", geo ch " << ix << " by " << dTshift[i];
1075 fhPullT_Smt_Off->SetBinContent(ix + 1, fhPullT_Smt_Off->GetBinContent(ix + 1) + dTshift[i] - dTMeanShift);
1076 }
1077 }
1078 else {
1079 iO[2] -= 2; // try other neighbor
1080 if (EvalDoublets(iLoc, iLoc + iO[1], iLoc + iO[2], dTshift)) { // returns vector of shifts
1081 double dTMeanShift = (dTshift[0] + dTshift[1] + dTshift[2]) / 3.;
1082 for (int i = 0; i < 3; i++) {
1083 int iStation = GetStationOfAddr(fDigiBdfPar->GetDetUId(iLoc + iO[i]));
1084 if (fiStationStatus[iStation] > 0) continue; // do not modify
1085 int ix = fMapRpcIdParInd[fDigiBdfPar->GetDetUId(iLoc + iO[i])]; // convert BDF to Geo counting
1086 LOG(info) << "UpdateDT- bdf ch " << iLoc + iO[i] << ", geo ch " << ix << " by " << dTshift[i];
1087 fhPullT_Smt_Off->SetBinContent(ix + 1,
1088 fhPullT_Smt_Off->GetBinContent(ix + 1) + dTshift[i] - dTMeanShift);
1089 }
1090 }
1091 iO[2] += 2; // restore offset
1092 }
1093 iO[2]--; // restore offset
1094 }
1095 }
1096 } break;
1097
1098 default: LOG(info) << "Correction mode not implemented!";
1099 }
1100
1101 if (NULL != fhPullT_Smt_Off) {
1102 // always extract residual widthes in T, X, Y, Z
1103 for (Int_t iStation = 0; iStation < static_cast<Int_t>(fMapRpcIdParInd.size()); iStation++) {
1104 TH1D* hResidualT = fhPullT_Smt->ProjectionY("_py", iStation + 1, iStation + 1);
1105 TH1D* hResidualX = fhPullX_Smt->ProjectionY("_py", iStation + 1, iStation + 1);
1106 TH1D* hResidualY = fhPullY_Smt->ProjectionY("_py", iStation + 1, iStation + 1);
1107 TH1D* hResidualZ = fhPullZ_Smt->ProjectionY("_py", iStation + 1, iStation + 1);
1108 if (hResidualT->GetEntries() > 100.) {
1109 Double_t dRMS = TMath::Abs(hResidualT->GetRMS());
1110 if (dRMS < RMSmin) dRMS = RMSmin;
1111 if (dRMS > 3. * fSIGT) dRMS = 3. * fSIGT;
1112 fhPullT_Smt_Width->SetBinContent(iStation + 1, dRMS);
1113 }
1114 if (hResidualX->GetEntries() > 100.) {
1115 Double_t dRMS = TMath::Abs(hResidualX->GetRMS());
1116 if (dRMS < 0.5 * fSIGX) dRMS = 0.5 * fSIGX;
1117 if (dRMS > 3. * fSIGX) dRMS = 3. * fSIGX;
1118 fhPullX_Smt_Width->SetBinContent(iStation + 1, dRMS);
1119 }
1120 if (hResidualY->GetEntries() > 100.) {
1121 Double_t dRMS = TMath::Abs(hResidualY->GetRMS());
1122 if (dRMS < 0.5 * fSIGY) dRMS = 0.5 * fSIGY;
1123 if (dRMS > 3. * fSIGY) dRMS = 3. * fSIGY;
1124 fhPullY_Smt_Width->SetBinContent(iStation + 1, dRMS);
1125 }
1126 if (hResidualZ->GetEntries() > 100.) {
1127 Double_t dRMS = TMath::Abs(hResidualZ->GetRMS());
1128 if (dRMS < 0.1) dRMS = 0.1;
1129 if (dRMS > 1.0) dRMS = 1.;
1130 fhPullZ_Smt_Width->SetBinContent(iStation + 1, dRMS);
1131 }
1132 }
1133
1134 fhPullT_Smt_Off->Write();
1135 fhPullX_Smt_Off->Write();
1136 fhPullY_Smt_Off->Write();
1137 fhPullZ_Smt_Off->Write();
1138 fhPullT_Smt_Width->Write();
1139 fhPullX_Smt_Width->Write();
1140 fhPullY_Smt_Width->Write();
1141 fhPullZ_Smt_Width->Write();
1142 }
1144 gFile = oldFile;
1145 gDirectory = oldDir;
1146
1147 fHist->Close();
1148
1149 return kTRUE;
1150}
1151
1152
1153// ----- Public method Exec --------------------------------------------
1154void CbmTofFindTracks::Exec(Option_t* opt)
1155{
1156 //if (fair::Logger::Logging(fair::Severity::debug)) { fDigiBdfPar->printParams(); }
1157 if (!fEventsColl) {
1158 // fTofHitArray = (TClonesArray*)fTofHitArrayIn->Clone();
1159 fTofHitArray = (TClonesArray*) fTofHitArrayIn;
1160 ExecFind(opt);
1161 }
1162 else {
1163 Int_t iNbTrks = 0;
1164 Int_t iNbCalHits = 0;
1165 fTrackArrayOut->Delete(); //Clear("C");
1166 fTofHitArrayOut->Delete(); //Clear("C");
1167 fTofUHitArrayOut->Delete(); //Clear("C");
1168 LOG(info) << "Process TS " << iTS << " with " << fEventsColl->GetEntriesFast() << " events";
1169 iTS++;
1170 for (Int_t iEvent = 0; iEvent < fEventsColl->GetEntriesFast(); iEvent++) {
1171 CbmEvent* tEvent = dynamic_cast<CbmEvent*>(fEventsColl->At(iEvent));
1172 LOG(debug) << "Process event " << iEvent << " with " << tEvent->GetNofData(ECbmDataType::kTofHit) << " hits";
1173 if ((tEvent->GetNofData(ECbmDataType::kTofHit)) < 1) continue;
1174
1175 if (fTofHitArray) fTofHitArray->Clear("C");
1176 Int_t iNbHits = 0;
1178 for (Int_t iHit = 0; iHit < (int) tEvent->GetNofData(ECbmDataType::kTofHit); iHit++) {
1179 Int_t iHitIndex = static_cast<Int_t>(tEvent->GetIndex(ECbmDataType::kTofHit, iHit));
1180 CbmTofHit* tHit = dynamic_cast<CbmTofHit*>(fTofHitArrayIn->At(iHitIndex));
1181 fTofHitIndexArray[iNbHits] = iHitIndex;
1182 new ((*fTofHitArray)[iNbHits++]) CbmTofHit(*tHit);
1183 }
1184
1185 ExecFind(opt, tEvent);
1186
1187 // --- In event-by-event mode: copy tracks to output array and register them to event
1188 LOG(info) << "Register " << fTrackArray->GetEntriesFast() << " tracks for event " << iEvent
1189 << ", starting at ind " << iNbTrks;
1190 for (Int_t iTrk = 0; iTrk < fTrackArray->GetEntriesFast(); iTrk++) {
1191 CbmTofTracklet* pTrk = (CbmTofTracklet*) fTrackArray->At(iTrk);
1192 new ((*fTrackArrayOut)[iNbTrks]) CbmTofTracklet(*pTrk);
1193 pTrk = (CbmTofTracklet*) fTrackArrayOut->At(iNbTrks);
1194 for (Int_t iHit = 0; iHit < pTrk->GetNofHits(); iHit++) //update to original index
1195 {
1196 pTrk->SetTofHitIndex(iHit, fTofHitIndexArray[pTrk->GetTofHitIndex(iHit)]);
1197 }
1198 tEvent->AddData(ECbmDataType::kTofTracklet, iNbTrks);
1199 iNbTrks++;
1200 }
1201 // Update TofHitArrayIn
1202 for (Int_t iHit = 0; iHit < fTofHitArray->GetEntriesFast(); iHit++) {
1203 CbmTofHit* tHit = dynamic_cast<CbmTofHit*>(fTofHitArray->At(iHit));
1204 new ((*fTofHitArrayOut)[iNbCalHits++]) CbmTofHit(*tHit);
1205 }
1206 // Add TofUHitArrayOut
1207 int iNbUHits = fTofUHitArrayOut->GetEntries();
1208 for (Int_t iHit = 0; iHit < fTofUHitArray->GetEntriesFast(); iHit++) {
1209 CbmTofHit* tHit = dynamic_cast<CbmTofHit*>(fTofUHitArray->At(iHit));
1210 new ((*fTofUHitArrayOut)[iNbUHits++]) CbmTofHit(*tHit);
1211 tEvent->AddData(ECbmDataType::kTofUHit, iNbUHits);
1212 }
1213
1214 fTrackArray->Delete();
1215 //fTrackArray->Clear();
1216 }
1217 }
1218}
1219
1220void CbmTofFindTracks::ExecFind(Option_t* /*opt*/, CbmEvent* tEvent)
1221{
1222 fiEvent++;
1224 if (NULL != fTofUHitArray) fTofUHitArray->Clear("C");
1225 if (NULL != fTrackArray) fTrackArray->Delete(); // reset
1226 //if (NULL != fTrackArray) fTrackArray->Clear(); // reset
1227 CbmTofHit* pBeamHit = NULL;
1228
1229 // recalibrate hits and count trackable hits
1230 for (Int_t iHit = 0; iHit < fTofHitArray->GetEntriesFast(); iHit++) {
1231 CbmTofHit* pHit = (CbmTofHit*) fTofHitArray->At(iHit);
1232 Int_t iDetId = (pHit->GetAddress() & DetMask);
1233
1234 Double_t dSIGX = GetSigX(iDetId);
1235 if (dSIGX == 0.) dSIGX = fSIGX;
1236 Double_t dSIGY = GetSigY(iDetId);
1237 if (dSIGY == 0.) dSIGY = fSIGY;
1238 Double_t dSIGZ = GetSigZ(iDetId);
1239 if (dSIGZ == 0.) dSIGZ = fSIGZ;
1240 TVector3 hitPosErr(dSIGX, dSIGY, dSIGZ); // include positioning uncertainty
1241 pHit->SetPositionError(hitPosErr);
1242
1243 Double_t dSIGT = GetSigT(iDetId);
1244 if (dSIGT == 0.) dSIGT = fSIGT;
1245 pHit->SetTimeError(dSIGT);
1246
1247 if (fiBeamCounter > -1) {
1248 // set diamond positions to (0,0,0) to allow inclusion in straight line fit
1249 if ((iDetId & DetMask) == fiBeamCounter) // modify diamond position
1250 {
1251 if (0. != fdBeamMomentumLab) {
1252 Double_t dTargetTimeOffset =
1253 pHit->GetZ() / fdBeamMomentumLab
1254 * TMath::Sqrt(TMath::Power(fdBeamMomentumLab, 2.) + TMath::Power(0.938271998, 2.)) / TMath::Ccgs() * 1.0e09;
1255 pHit->SetTime(pHit->GetTime() - dTargetTimeOffset);
1256 }
1257
1258 TVector3 hitPos(0., 0., 0.);
1259 // TVector3 hitPos(pHit->GetX(), pHit->GetY(), 0.);
1260 // TVector3 hitPosErr(1.,1.,5.0); // including positioning uncertainty
1261 pHit->SetPosition(hitPos);
1262 TVector3 hitPosErr0(1., 1., 1.); // including positioning uncertainty
1263 pHit->SetPositionError(hitPosErr0); //
1264 if (NULL != pBeamHit) {
1265 if (pHit->GetTime() < pBeamHit->GetTime()) pBeamHit = pHit;
1266 }
1267 else
1268 pBeamHit = pHit;
1269 }
1270 }
1271
1273 Int_t iHitAddress = pHit->GetAddress();
1274 Int_t iModuleType = CbmTofAddress::GetSmType(iHitAddress);
1275 Int_t iModuleIndex = CbmTofAddress::GetSmId(iHitAddress);
1276 Int_t iCounterIndex = CbmTofAddress::GetRpcId(iHitAddress);
1277
1278 CbmTofCell* fChannelInfo = fDigiPar->GetCell(iHitAddress);
1279
1280 Double_t dSignalPropagationTime =
1281 0.5
1282 * (fChannelInfo->GetSizey() >= fChannelInfo->GetSizex() ? fChannelInfo->GetSizey() : fChannelInfo->GetSizex())
1283 / fDigiBdfPar->GetSigVel(iModuleType, iModuleIndex, iCounterIndex);
1284
1285 pHit->SetTime(pHit->GetTime() - dSignalPropagationTime);
1286 }
1287
1288 // tune positions and times
1289 if ((iDetId & DetMask) != fiBeamCounter) { // do not modify diamond position
1290 Int_t iRpcInd = fMapRpcIdParInd[iDetId];
1291 pHit->SetTime(pHit->GetTime() + fvToff[iRpcInd]);
1292 pHit->SetX(pHit->GetX() + fvXoff[iRpcInd]);
1293 pHit->SetY(pHit->GetY() + fvYoff[iRpcInd]);
1294 pHit->SetZ(pHit->GetZ() + fvZoff[iRpcInd]);
1295 }
1296
1297 Int_t iSt = GetStationOfAddr(iDetId);
1298 if (iSt >= (Int_t) fStationHMul.size()) {
1299 LOG(fatal) << Form("Invalid station # %d for detId 0x%08x, TSR %d%d%d", iSt, iDetId,
1301 CbmTofAddress::GetRpcId(iDetId));
1302 }
1303
1304 if ((Int_t) fMapStationRpcId.size() > fNTofStations) {
1305 PrintSetup();
1306 LOG(fatal) << "Invalid NTofStations " << fNTofStations << ", " << fMapStationRpcId.size();
1307 }
1308
1309 LOG(debug) << Form("ExecFind Hit %02d, addr 0x%08x, TSR %d%d%d, sta %02d, %02d, HM "
1310 "%02d, X %6.2f(%3.2f) Y "
1311 "%6.2f(%3.2f) Z %6.2f(%3.2f) T %12.2f(%3.2f)",
1312 iHit, pHit->GetAddress(), CbmTofAddress::GetSmType(iDetId), CbmTofAddress::GetSmId(iDetId),
1314 fStationHMul[GetStationOfAddr(iDetId)], pHit->GetX(), pHit->GetDx(), pHit->GetY(), pHit->GetDy(),
1315 pHit->GetZ(), pHit->GetDz(), pHit->GetTime(), pHit->GetTimeError());
1316
1317 MarkStationFired(iSt);
1318 }
1319
1320 /* call moved to Clusterizer
1321 if (fiCalOpt > 99) {
1322 //LOG(info) << "Call TofCalibrator with pBeam "<<pBeamHit;
1323 if (NULL != pBeamHit) fTofCalibrator->FillHitCalHist(pBeamHit, fiCalOpt, tEvent, fTofHitArrayIn);
1324 }
1325 */
1326
1327 LOG(debug) << Form("CbmTofFindTracks::Exec NStationsFired %d >= %d Min ?, NbStations %d", GetNStationsFired(),
1329
1330 if (GetNStationsFired() < GetMinNofHits()) {
1331 fInspectEvent = kFALSE; // mark event as non trackable
1332 }
1333 else
1334 fInspectEvent = kTRUE;
1335
1336 CheckMaxHMul();
1337 // resort Hit array with respect to time, FIXME danger: links to digis become invalid (???, check!!!)
1338 // fTofHitArray->Sort(fTofHitArray->GetEntriesFast()); // feature not available
1339
1340 if (fInspectEvent && fNTofStations > 1) {
1341 fStart.Set();
1342 //fTrackArray->Clear("C+C");
1344 // fTrackArray->Compress();
1345 fStop.Set();
1346 fdTrackingTime = fStop.GetSec() - fStart.GetSec() + (fStop.GetNanoSec() - fStart.GetNanoSec()) / 1e9;
1347
1348 LOG(debug) << Form("CbmTofFindTracks::Exec found %d Tracklets in %f sec", fTrackArray->GetEntriesFast(),
1350
1351 FindVertex();
1352
1353 if (fbDoHistos) FillHistograms(tEvent);
1354 }
1355
1356 FillUHits(); // put unused hits into TClonesArray
1357}
1358// -------------------------------------------------------------------------
1359
1360
1361// ----- Public method Finish ------------------------------------------
1363{
1364 if (fiEvent < 2) return; // preserve calibration histos in event display
1365 FairLogger::GetLogger()->SetLogScreenLevel("info");
1367 if (fbDoHistos) WriteHistos();
1368
1369 LOG(info) << Form(" CbmTofFindTracks::Finished for opt %d", fiCalOpt);
1370}
1371// -------------------------------------------------------------------------
1372
1374{
1375
1376 TDirectory* oldir = gDirectory; // <= To prevent histos from being sucked in by the param file of the TRootManager!
1377 gROOT->cd(); // <= To prevent histos from being sucked in by the param file of the TRootManager !
1378
1379 // define histos here
1380
1381 Double_t nSmt = fMapRpcIdParInd.size();
1382 LOG(info) << Form(" CbmTofFindTracks::CreateHistograms for %d counters, %d stations ", (Int_t) nSmt, fNTofStations);
1383
1384 fhTrklMul = new TH1F(Form("hTrklMul"), Form("Tracklet Multiplicity; MulTracklet"), 100, 0, 100);
1385
1387 new TH1F(Form("hAllHitsStation"), Form("Reconstructed Hits; Station #"), fNTofStations, 0, fNTofStations);
1388 fhAllHitsSmTypes = new TH1F(Form("hAllHitsSmTypes"), Form("Reconstructed Hits; SmType #"), 10, 0, 10);
1389
1390 fhUsedHitsStation = new TH1F(Form("hUsedHitsStation"), Form("Used (HMul>2) / Reconstructed Hits; Station #"),
1392
1393 fhTrklChi2 = new TH2F(Form("hTrklChi2"), Form("Tracklet Chi; HMul_{Tracklet}; #chi"), 8, 2, 10, 100, 0,
1394 ((CbmTofTrackFinderNN*) fFinder)->GetChiMaxAccept());
1395
1397 new TH2F(Form("hTrackingTimeNhits"), Form("Tracking Time; NHits; #Deltat (s)"), 200, 0, 200, 50, 0, 0.2);
1398
1400 new TH2F(Form("hTrklMulNhits"), Form("Tracklet Multiplicity; NHits; NTracklet"), 200, 0, 200, 30, 0, 30);
1401
1403 new TH2F(Form("hTrklMulMaxMax-1"), Form("Tracklet Multiplicity; TMulMax; TMulMax-1"), 10, 0, 10, 10, 0, 10);
1404 fhTrklMul3D =
1405 new TH3F(Form("hTrklMul3D"), Form("Tracklet Multiplicities; TMul3; TMul4; TMul5"), 10, 0, 10, 10, 0, 10, 10, 0, 10);
1406
1407 fhTrklHMul =
1408 new TH2F(Form("hTrklHMul"), Form("Tracklet Hit - Multiplicity; HMul_{Tracklet}; Mul_{HMul}"), 8, 2, 10, 20, 0, 20);
1409 fhTrklZ0xHMul = new TH2F(Form("hTrklZ0xHMul"), Form("Tracklet Z0x vs. Hit - Multiplicity; HMul_{Tracklet}; Z0x"), 8,
1410 2, 10, 100, -500, 500);
1411 fhTrklZ0yHMul = new TH2F(Form("hTrklZ0yHMul"), Form("Tracklet Z0y vs. Hit - Multiplicity; HMul_{Tracklet}; Z0y"), 8,
1412 2, 10, 100, -300, 300);
1413
1414 double dTxLim = ((CbmTofTrackFinderNN*) fFinder)->GetTxLIM();
1415 double dTyLim = ((CbmTofTrackFinderNN*) fFinder)->GetTyLIM();
1416 double dTxMean = ((CbmTofTrackFinderNN*) fFinder)->GetTxMean();
1417 double dTyMean = ((CbmTofTrackFinderNN*) fFinder)->GetTyMean();
1418 double dTxMin = dTxMean - dTxLim;
1419 double dTxMax = dTxMean + dTxLim;
1420 double dTyMin = dTyMean - dTyLim;
1421 double dTyMax = dTyMean + dTyLim;
1422
1423 fhTrklTxHMul = new TH2F(Form("hTrklTxHMul"), Form("Tracklet Tx vs. Hit - Multiplicity; HMul_{Tracklet}; Tx"), 8, 2,
1424 10, 100, dTxMin, dTxMax);
1425 fhTrklTyHMul = new TH2F(Form("hTrklTyHMul"), Form("Tracklet Ty vs. Hit - Multiplicity; HMul_{Tracklet}; Ty"), 8, 2,
1426 10, 100, dTyMin, dTyMax);
1427
1428 fhTrklTyTx = new TH2F(Form("hTrklTyTx"), Form("Tracklet Ty vs. Tx; Tx; Ty"), 50, dTxMin, dTxMax, 50, dTyMin, dTyMax);
1429
1430 Double_t TTMAX = 0.2;
1431 fhTrklTtHMul = new TH2F(Form("hTrklTtHMul"), Form("Tracklet Tt vs. Hit - Multiplicity; HMul_{Tracklet}; Tt"), 8, 2,
1432 10, 100, -TTMAX, TTMAX);
1434 new TH2F(Form("hTrklVelHMul"), Form("Tracklet Vel vs. Hit - Multiplicity; HMul_{Tracklet}; v (cm/ns)"), 8, 2, 10,
1435 100, 0., 50.);
1436 fhTrklT0HMul = new TH2F(Form("hTrklT0HMul"), Form("Tracklet T0 vs. Hit - Multiplicity; HMul_{Tracklet}; T0"), 8, 2,
1437 10, 100, -0.5, 0.5);
1438
1439 fhTrklT0Mul =
1440 new TH2F(Form("hTrklT0Mul"), Form("Tracklet #DeltaT0 vs. Trkl - Multiplicity; Mul_{Tracklet}; #Delta(T0)"), 10, 0,
1441 10, 100, -2., 2.);
1442 fhTrklDT0SmMis = new TH2F(Form("hTrklDT0SmMis"), Form("Tracklet DeltaT0 vs. Trkl - ID; SmType_{missed}; #Delta(T0)"),
1443 10, 0, 10, 100, -2., 2.);
1445 new TH2F(Form("hTrklDT0SmMis2"), Form("Tracklet DeltaT0 vs. Station - ID; St2_{missed}; #Delta(T0)"), 50, 0, 50,
1446 100, -2., 2.);
1447
1448 Double_t X0MAX = 40.;
1449 fhTrklXY0_0 = new TH2F(Form("hTrklXY0_0"), Form("Tracklet XY at z=0 for hmulmax ; x (cm); y (cm)"), 100, -X0MAX,
1450 X0MAX, 100, -X0MAX, X0MAX);
1451 fhTrklXY0_1 = new TH2F(Form("hTrklXY0_1"), Form("Tracklet XY at z=0 for hmulmax-1 ; x (cm); y (cm)"), 100, -X0MAX,
1452 X0MAX, 100, -X0MAX, X0MAX);
1453 fhTrklXY0_2 = new TH2F(Form("hTrklXY0_2"), Form("Tracklet XY at z=0 for hmulmax-2 ; x (cm); y (cm)"), 100, -X0MAX,
1454 X0MAX, 100, -X0MAX, X0MAX);
1455
1456 fhTrklX0_TX = new TH2F(Form("hTrklX0_TX"), Form("Tracklet X0 vs TX at z=0 for hmulmax ; tx (); x0 (cm)"), 100, -0.4,
1457 0.4, 100, -X0MAX, X0MAX);
1458 fhTrklY0_TX = new TH2F(Form("hTrklY0_TX"), Form("Tracklet Y0 vs TX at z=0 for hmulmax ; tx (); y0 (cm)"), 100, -0.4,
1459 0.4, 100, -X0MAX, X0MAX);
1460 fhTrklX0_TY = new TH2F(Form("hTrklX0_TY"), Form("Tracklet X0 vs TY at z=0 for hmulmax ; ty (); x0 (cm)"), 100, -0.4,
1461 0.4, 100, -X0MAX, X0MAX);
1462 fhTrklY0_TY = new TH2F(Form("hTrklY0_TY"), Form("Tracklet Y0 vs TY at z=0 for hmulmax ; ty (); y0 (cm)"), 100, -0.4,
1463 0.4, 100, -X0MAX, X0MAX);
1464
1465 Double_t DT0MAX = 5.;
1466 if (fT0MAX == 0) fT0MAX = DT0MAX;
1467 fhPullT_Smt = new TH2F(Form("hPullT_Smt"), Form("Tracklet ResiT vs RpcInd ; RpcInd ; #DeltaT (ns)"), nSmt, 0, nSmt,
1468 501, -fT0MAX, fT0MAX);
1469 Double_t DX0MAX = 5.;
1470 fhPullX_Smt = new TH2F(Form("hPullX_Smt"), Form("Tracklet ResiX vs RpcInd ; RpcInd ; #DeltaX (cm)"), nSmt, 0, nSmt,
1471 100, -DX0MAX, DX0MAX);
1472 Double_t DY0MAX = 5.;
1473 fhPullY_Smt = new TH2F(Form("hPullY_Smt"), Form("Tracklet ResiY vs RpcInd ; RpcInd ; #DeltaY (cm)"), nSmt, 0, nSmt,
1474 100, -DY0MAX, DY0MAX);
1475 Double_t DZ0MAX = 20.;
1476 fhPullZ_Smt = new TH2F(Form("hPullZ_Smt"), Form("Tracklet ResiZ vs RpcInd ; RpcInd ; #DeltaZ (cm)"), nSmt, 0, nSmt,
1477 100, -DZ0MAX, DZ0MAX);
1478
1479 fhTOff_Smt = new TH2F(Form("hTOff_Smt"), Form("Tracklet TOff; RpcInd ; #DeltaTOff (ns)"), nSmt, 0, nSmt, 501,
1480 -fT0MAX * 5, fT0MAX * 5);
1481 fhTOff_HMul2 = new TH2F(Form("hTOff_HMul2"), Form("Tracklet TOff(HMul2); RpcInd ; TOff (ns)"), nSmt, 0, nSmt, 500,
1482 -fT0MAX, fT0MAX);
1483
1484 Double_t DTTMAX = 0.11;
1485 fhDeltaTt_Smt = new TH2F(Form("hDeltaTt_Smt"), Form("Tracklet DeltaTt; RpcInd ; #DeltaTt (ns/cm)"), nSmt, 0, nSmt,
1486 100, -DTTMAX, DTTMAX);
1487
1489 new TH2F(Form("hDeltaTc_Smt"), Form("Tracklet DeltaTc; RpcInd ; #DeltaTc (ns)"), nSmt, 0, nSmt, 160, -2., 2.);
1490
1491 vhPullX.resize(fNTofStations);
1492 vhPullY.resize(fNTofStations);
1493 vhPullZ.resize(fNTofStations);
1494 vhPullT.resize(fNTofStations);
1495 vhPullTB.resize(fNTofStations);
1496 vhTrefRms.resize(fNTofStations);
1497 vhFitDT0.resize(fNTofStations);
1498 vhFitT0Err.resize(fNTofStations);
1499 vhFitTt.resize(fNTofStations);
1500 vhFitTtErr.resize(fNTofStations);
1501 vhFitDTMean.resize(fNTofStations);
1509 vhXY_DX.resize(fNTofStations);
1510 vhXY_DY.resize(fNTofStations);
1511 vhXY_DT.resize(fNTofStations);
1512 vhXY_TOT.resize(fNTofStations);
1513 vhXY_CSZ.resize(fNTofStations);
1514 vhUDXDY_DT.resize(fNTofStations);
1515 vhUCDXDY_DT.resize(fNTofStations);
1516
1517 for (Int_t iSt = 0; iSt < fNTofStations; iSt++) {
1518 vhPullX[iSt] =
1519 new TH1F(Form("hPullX_Station_%d", iSt), Form("hResiX_Station_%d; #DeltaX (cm)", iSt), 99, -DX0MAX, DX0MAX);
1520 vhPullY[iSt] =
1521 new TH1F(Form("hPullY_Station_%d", iSt), Form("hResiY_Station_%d; #DeltaY (cm)", iSt), 99, -DY0MAX, DY0MAX);
1522 vhPullZ[iSt] =
1523 new TH1F(Form("hPullZ_Station_%d", iSt), Form("hResiZ_Station_%d; #DeltaZ (cm)", iSt), 99, -50., 50.);
1524 vhPullT[iSt] =
1525 new TH1F(Form("hPullT_Station_%d", iSt), Form("hResiT_Station_%d; #DeltaT (ns)", iSt), 59, -fT0MAX, fT0MAX);
1526 vhPullTB[iSt] = new TH1F(Form("hPullTB_Station_%d", iSt), Form("hResiTB_Station_%d; #DeltaT (ns)", iSt), 59,
1527 -1.25 * fT0MAX, 1.25 * fT0MAX);
1528 vhTrefRms[iSt] = new TH1F(Form("hTrefRms_Station_%d", iSt), Form("hTrefRms_Station_%d; RMS (ns)", iSt), 40, 0, 2.);
1529 vhFitDT0[iSt] =
1530 new TH1F(Form("hFitDT0_Station_%d", iSt), Form("hFitDT0_Station_%d; #Deltat (ns)", iSt), 40, -5, 5.);
1531 vhFitT0Err[iSt] =
1532 new TH1F(Form("hFitT0Err_Station_%d", iSt), Form("hFitT0Err_Station_%d; rel. error ()", iSt), 40, 0, 0.1);
1533 vhFitTt[iSt] =
1534 new TH1F(Form("hFitTt_Station_%d", iSt), Form("hFitTt_Station_%d; inverse velocity (ns/cm)", iSt), 40, 0, 0.1);
1535 vhFitTtErr[iSt] = new TH1F(Form("hFitTtErr_Station_%d", iSt),
1536 Form("hFitTtErr_Station_%d; inverse velocity error()", iSt), 40, 0, 0.4);
1537 vhFitDTMean[iSt] = new TH1F(Form("hFitDTMean_Station_%d", iSt),
1538 Form("hFitDTMean_Station_%d; extrapolation distance (ns)", iSt), 40, -5, 5.);
1539 vhFitDTMeanErr[iSt] = new TH1F(Form("hFitDTMeanErr_Station_%d", iSt),
1540 Form("hFitDTMeanErr_Station_%d; extrapolation error (ns)", iSt), 40, 0, 0.4);
1541
1542 const Double_t TOTmax = 50.;
1543 vhResidualTBWalk[iSt] =
1544 new TH2F(Form("hResidualTBWalk_Station_%d", iSt), Form("hResidualTBWalk_Station_%d; #DeltaT (ns)", iSt), TOTmax,
1545 0., TOTmax, 59, -1.25 * fT0MAX, 1.25 * fT0MAX);
1546 vhResidualYWalk[iSt] =
1547 new TH2F(Form("hResidualYWalk_Station_%d", iSt), Form("hResidualYWalk_Station_%d; #DeltaT (ns)", iSt), TOTmax,
1548 0., TOTmax, 59, -DY0MAX, DY0MAX);
1549 Double_t XSIZ = 16.;
1550 Int_t Nbins = 32.;
1551 CbmTofCell* fChannelInfo = fDigiPar->GetCell(fMapStationRpcId[iSt]);
1552 if (NULL == fChannelInfo) {
1553 LOG(fatal) << "Geometry for station " << iSt << Form(", Rpc 0x%08x ", fMapStationRpcId[iSt]) << " not defined ";
1554 return;
1555 }
1556 Int_t NbinsX = 2 * XSIZ / fChannelInfo->GetSizex();
1557 vhXY_AllTracks[iSt] = new TH2F(Form("hXY_AllTracks_%d", iSt), Form("hXY_AllTracks_%d; x(cm); y (cm)", iSt), NbinsX,
1558 -XSIZ, XSIZ, Nbins, -XSIZ, XSIZ);
1559 vhXY_AllStations[iSt] = new TH2F(Form("hXY_AllStations_%d", iSt), Form("hXY_AllStations_%d; x(cm); y (cm)", iSt),
1560 NbinsX, -XSIZ, XSIZ, Nbins, -XSIZ, XSIZ);
1561 vhXY_AllFitStations[iSt] =
1562 new TH2F(Form("hXY_AllFitStations_%d", iSt), Form("hXY_AllFitStations_%d; x(cm); y (cm)", iSt), NbinsX, -XSIZ,
1563 XSIZ, Nbins, -XSIZ, XSIZ);
1564 vhXY_MissedStation[iSt] =
1565 new TH2F(Form("hXY_MissedStation_%d", iSt), Form("hXY_MissedStation_%d; x(cm); y (cm)", iSt), NbinsX, -XSIZ,
1566 XSIZ, Nbins, -XSIZ, XSIZ);
1567 vhXY_DX[iSt] = new TH3F(Form("hXY_DX_%d", iSt), Form("hXY_DX_%d; x(cm); y (cm); #DeltaX (cm)", iSt), NbinsX, -XSIZ,
1568 XSIZ, Nbins, -XSIZ, XSIZ, Nbins, -2., 2.);
1569 vhXY_DY[iSt] = new TH3F(Form("hXY_DY_%d", iSt), Form("hXY_DY_%d; x(cm); y (cm); #DeltaY (cm)", iSt), NbinsX, -XSIZ,
1570 XSIZ, Nbins, -XSIZ, XSIZ, Nbins, -2., 2.);
1571 vhXY_DT[iSt] = new TH3F(Form("hXY_DT_%d", iSt), Form("hXY_DT_%d; x(cm); y (cm); #DeltaT (ns)", iSt), NbinsX, -XSIZ,
1572 XSIZ, Nbins, -XSIZ, XSIZ, Nbins, -0.5, 0.5);
1573 vhXY_TOT[iSt] = new TH3F(Form("hXY_TOT_%d", iSt), Form("hXY_TOT_%d; x(cm); y (cm); TOT (a.u.)", iSt), NbinsX,
1574 -XSIZ, XSIZ, Nbins, -XSIZ, XSIZ, Nbins, 0., 20.);
1575 vhXY_CSZ[iSt] = new TH3F(Form("hXY_CSZ_%d", iSt), Form("hXY_CSZ_%d; x(cm); y (cm); CSZ ()", iSt), NbinsX, -XSIZ,
1576 XSIZ, Nbins, -XSIZ, XSIZ, 6, 1., 7.);
1577 vhUDXDY_DT[iSt] = new TH3F(Form("hUDXDY_DT_%d", iSt),
1578 Form("Unused missing hit - DXDY_DT_%d; #Deltax "
1579 "(cm); #Deltay (cm); #DeltaT (ns)",
1580 iSt),
1581 11, -3., 3., 11, -3., 3., 101, -50., 50.);
1582 vhUCDXDY_DT[iSt] = new TH3F(Form("hUCDXDY_DT_%d", iSt),
1583 Form("Unused close hit - DXDY_DT_%d; #Deltax "
1584 "(cm); #Deltay (cm); #DeltaT (ns)",
1585 iSt),
1586 11, -3., 3., 11, -3., 3., 101, -50., 50.);
1587 }
1588
1589
1590 // vertex histograms
1591 Double_t NNORM = 40.;
1592 fhVTXNorm = new TH1F(Form("hVTXNorm"), Form("Vertex Normalisation; #_{TrackletHits}"), NNORM, 0, NNORM);
1593 fhVTX_XY0 =
1594 new TH2F(Form("hVTX_XY0"), Form("Vertex XY at z=0 ; x (xm); y (cm)"), 100, -X0MAX, X0MAX, 100, -X0MAX, X0MAX);
1595 fhVTX_DT0_Norm = new TH2F(Form("hVTX_DT0_Norm"), Form("Vertex #DeltaT at z=0 ; #_{TrackletHits}; #DeltaT (ns)"),
1596 NNORM, 0, NNORM, 100, -2., 2.);
1597
1598 gDirectory->cd(oldir->GetPath()); // <= To prevent histos from being sucked in by the param file of the TRootManager!
1599}
1600
1602{
1603 fVTX_T = 0.; //reset
1604 fVTX_X = 0.;
1605 fVTX_Y = 0.;
1606 fVTX_Z = 0.;
1607 fVTXNorm = 0.;
1608
1609 for (Int_t iTrk = 0; iTrk < fTrackArray->GetEntriesFast(); iTrk++) {
1610 CbmTofTracklet* pTrk = (CbmTofTracklet*) fTrackArray->At(iTrk);
1611 if (NULL == pTrk) continue;
1612 Double_t w = pTrk->GetNofHits();
1613 LOG(debug1) << Form("CbmTofFindTracks::FindVertex: N %3.0f, w %3.0f, min %d", fVTXNorm, w, fMinNofHits);
1614
1615 if (w > (Double_t) fMinNofHits) { // for further analysis request minimum number of hits
1616 fVTXNorm += w;
1617 fVTX_T += w * pTrk->GetFitT(0.);
1618 fVTX_X += w * pTrk->GetFitX(0.);
1619 fVTX_Y += w * pTrk->GetFitY(0.);
1620 }
1621 }
1622 if (fVTXNorm > 0.) {
1623 fVTX_T /= fVTXNorm;
1624 fVTX_X /= fVTXNorm;
1625 fVTX_Y /= fVTXNorm;
1626 fVTX_Z /= fVTXNorm;
1627 }
1628 LOG(debug) << Form("CbmTofFindTracks::FindVertex: N %3.0f, T %6.2f, X=%6.2f, Y=%6.2f Z=%6.2f ", fVTXNorm, fVTX_T,
1629 fVTX_X, fVTX_Y, fVTX_Z);
1630}
1631
1632static Int_t iWarnNotDefined = 0;
1633
1635{
1636 // Locate reference ("beam counter") hit
1637 CbmTofHit* pRefHit = NULL;
1638 Double_t RefMinTime = 1.E300;
1639 for (Int_t iHit = 0; iHit < fTofHitArray->GetEntriesFast(); iHit++) { // loop over Hits
1640 CbmTofHit* pHit = (CbmTofHit*) fTofHitArray->At(iHit);
1641 Int_t iAddr = (pHit->GetAddress() & DetMask);
1642 if (fiBeamCounter != -1) {
1643 if (iAddr == fiBeamCounter)
1644 if (pHit->GetTime() < RefMinTime) {
1645 pRefHit = pHit;
1646 RefMinTime = pRefHit->GetTime();
1647 }
1648 }
1649 else { // take earliest hit as reference
1650 if (pHit->GetTime() < RefMinTime) {
1651 pRefHit = pHit;
1652 RefMinTime = pRefHit->GetTime();
1653 }
1654 }
1655 }
1656 if (fiBeamCounter != -1 && NULL == pRefHit) return;
1657 //LOG(debug)<<"FillH_start";
1658 //PrintSetup();
1659
1660 std::vector<Int_t> HMul;
1661 HMul.resize(fNTofStations + 1);
1662 // HMul.clear();
1663
1664 fhTrklMul->Fill(fTrackArray->GetEntriesFast());
1665
1666 Int_t iTMul = 0;
1667 for (Int_t iTrk = 0; iTrk < fTrackArray->GetEntriesFast(); iTrk++) {
1668 CbmTofTracklet* pTrk = (CbmTofTracklet*) fTrackArray->At(iTrk);
1669 if (NULL == pTrk) continue;
1670 if (pTrk->GetNofHits() > fNTofStations + 1) { // +1 for virtual T0
1671 LOG(error) << "CbmTofFindTracks::FillHistograms: more hits (" << pTrk->GetNofHits() << ") than stations ("
1672 << fNTofStations << ")";
1673 continue;
1674 }
1675
1676 //if (pTrk->GetTrackTy()>0.) continue; // Debugging!!!!
1677 //LOG(debug)<<"FillH_HMul "<<iTrk<<", "<<pTrk->GetNofHits();
1678 //PrintSetup();
1679
1680 HMul[pTrk->GetNofHits()]++;
1681
1682 if (pTrk->GetNofHits() >= 2) { // initial offset calibration
1683 //Int_t iH0 = pTrk->GetStationHitIndex(fMapStationRpcId[0]); // Hit index for station 0 (Diamond)
1684 //if(iH0<0) continue; // Station 0 not part of tracklet
1685 // Int_t iDetId0 = pTrk->GetTofDetIndex(0); // DetId of 1. Hit (FU) not used
1686 // Int_t iSt0 = GetStationOfAddr(iDetId0); // Station of 1. Hit (FU) not used
1687 CbmTofHit* pHit0 = pTrk->GetTofHitPointer(0);
1688 Double_t dTRef0 = pHit0->GetTime() - pHit0->GetR() * fTtTarg;
1689
1690 for (Int_t iH = 1; iH < pTrk->GetNofHits(); iH++) {
1691 Int_t iDetId = pTrk->GetTofDetIndex(iH); // DetId of iH. Hit
1692 CbmTofHit* pHit = pTrk->GetTofHitPointer(iH);
1693 Int_t iSt = GetStationOfAddr(iDetId); // Station of 1. Hit
1694 if (iSt < 0) continue;
1695 Double_t dTOff = pHit->GetTime() - pHit->GetR() * fTtTarg - dTRef0;
1696 LOG(debug1) << Form("<D> CbmTofFindTracks::FillHistograms: iDetId1 "
1697 "0x%08x, iST1 = %d with dTOff %f at RpcInd %d",
1698 iDetId, iSt, dTOff, fMapRpcIdParInd[fMapStationRpcId[iSt]]);
1699 fhTOff_HMul2->Fill((Double_t) fMapRpcIdParInd[fMapStationRpcId[iSt]], dTOff);
1700 } // loop over tracklets' hits
1701 }
1702
1703 if (pTrk->GetNofHits() >= fMinNofHits) { // for further analysis request at least 3 matched hits
1704 iTMul++;
1705 fhTrklChi2->Fill(pTrk->GetNofHits(), pTrk->GetChiSq());
1706
1707 fhTrklChi2->Fill(pTrk->GetNofHits(), pTrk->GetChiSq());
1708
1709 if (fiCalOpt > 0) // && fiCalOpt < 10 )
1710 if (pTrk->GetTt() > fdTtMin) fTofCalibrator->FillCalHist(pTrk, fiCalOpt, tEvent);
1711
1712 CbmTofTrackletParam* tPar = pTrk->GetTrackParameter();
1713 Double_t dTt = pTrk->GetTt();
1714 LOG(debug) << Form("Trk %d info: Lz=%6.2f Z0x=%6.2f Z0y=%6.2f Tt=%6.4f, ", iTrk, tPar->GetLz(), pTrk->GetZ0x(),
1715 pTrk->GetZ0y(), dTt)
1716 << tPar->ToString();
1717
1718 fhTrklZ0xHMul->Fill(pTrk->GetNofHits(), pTrk->GetFitX(0.));
1719 fhTrklZ0yHMul->Fill(pTrk->GetNofHits(), pTrk->GetFitY(0.));
1720 fhTrklTxHMul->Fill(pTrk->GetNofHits(), tPar->GetTx());
1721 fhTrklTyHMul->Fill(pTrk->GetNofHits(), tPar->GetTy());
1722 fhTrklTyTx->Fill(tPar->GetTx(), tPar->GetTy());
1723 fhTrklTtHMul->Fill(pTrk->GetNofHits(), dTt);
1724
1725 if (kFALSE) { //kFALSE && fdTtMin < -10.) { // print event number for inspection in event display
1726 //Double_t zPosTar = 160.;
1727 //if (TMath::Abs(tPar->GetTx())<0.2 && TMath::Abs(tPar->GetTy())<0.1) {
1728 //if(TMath::Abs(pTrk->GetFitX(0.))<5 && TMath::Abs(pTrk->GetFitY(0.))<5) {
1729 //if ((TMath::Abs(pTrk->GetFitX(zPosTar) - 20.) < 10 && TMath::Abs(pTrk->GetFitY(zPosTar) - 15.) < 5)
1730 if (pTrk->GetNofHits() > 8) {
1731 LOG(info) << Form("Found tracklet candidate (x0,y0)=(%5.2f,%5.2f), Hmul %d, v %6.2f in event %d ",
1732 pTrk->GetFitX(0.), pTrk->GetFitY(0.), pTrk->GetNofHits(), 1. / dTt, fiEvent)
1733 << " within " << fTofHitArray->GetEntriesFast() << " hits ";
1734 }
1735 }
1736 switch (GetNReqStations() - pTrk->GetNofHits()) {
1737 case 0: // max hit number
1738 fhTrklXY0_0->Fill(pTrk->GetFitX(0.), pTrk->GetFitY(0.));
1739 fhTrklX0_TX->Fill(tPar->GetTx(), pTrk->GetFitX(0.));
1740 fhTrklY0_TX->Fill(tPar->GetTx(), pTrk->GetFitY(0.));
1741 fhTrklX0_TY->Fill(tPar->GetTy(), pTrk->GetFitX(0.));
1742 fhTrklY0_TY->Fill(tPar->GetTy(), pTrk->GetFitY(0.));
1743 break;
1744 case 1: fhTrklXY0_1->Fill(pTrk->GetFitX(0.), pTrk->GetFitY(0.)); break;
1745 case 2: fhTrklXY0_2->Fill(pTrk->GetFitX(0.), pTrk->GetFitY(0.)); break;
1746 default:;
1747 }
1748
1749 if (fiBeamCounter > -1 && fdR0Lim > 0.) // consider only tracks originating from nominal interaction point
1750 {
1751 if (pTrk->GetR0() > fdR0Lim) continue;
1752 }
1753
1754 if (dTt > fdTtMin)
1755 for (Int_t iSt = 0; iSt < fNTofStations; iSt++) {
1756 Int_t iH = pTrk->GetStationHitIndex(fMapStationRpcId[iSt]); // Station Hit index
1757 if (iH < 0) continue; // Station not part of tracklet
1758 fhUsedHitsStation->Fill(iSt);
1759
1760 if (pTrk->GetNofHits() < GetNReqStations()) continue; // fill Pull histos only for complete tracks
1761 CbmTofHit* pHit = (CbmTofHit*) fTofHitArray->At(iH);
1762 //if (0 == fMapStationRpcId[iSt]) pHit->SetTime(pTrk->GetT0()); // set time of fake hit, abandoned
1763 /*
1764 LOG(info) << " -D- CbmTofFindTracks::FillHistograms: "<< iSt <<", "
1765 <<fMapStationRpcId[iSt]<<", "<< iH <<", "<< iH0 <<", "<<pHit->ToString();
1766 */
1767 Double_t dDZ = pHit->GetZ() - tPar->GetZ(); // z- Distance to reference point
1768 Double_t dDX = pHit->GetX() - pTrk->GetFitX(pHit->GetZ()); // - tPar->GetX() - tPar->GetTx()*dDZ;
1769 Double_t dDY = pHit->GetY() - pTrk->GetFitY(pHit->GetZ()); // - tPar->GetTy()*dDZ;
1770 Double_t dDT = pHit->GetTime() - pTrk->GetFitT(pHit->GetZ()); // pTrk->GetTdif(fMapStationRpcId[iSt]);
1771 Double_t dTexp = fTrackletTools->GetTexpected(pTrk, fMapStationRpcId[iSt], pHit, fTtLight);
1772 // fTrackletTools->GetTdif(pTrk, fMapStationRpcId[iSt],pHit); // ignore pHit in calc of reference
1773 Double_t dDTB = pHit->GetTime() - dTexp;
1774 Double_t dTexpErr = fTrackletTools->GetTexpectedError(pTrk, fMapStationRpcId[iSt], pHit, dTexp);
1775 Double_t dTOT = pHit->GetCh() / 10.; // misuse of channel field
1776
1777 Double_t dZZ = pHit->GetZ() - tPar->GetZy(pHit->GetY());
1778 LOG(debug) << Form(" St %d Id 0x%08x Hit %2d, Z %6.2f - DX %6.2f, DY %6.2f, "
1779 "Z %6.2f, DT %6.2f, %6.2f, ZZ %6.2f, Tt %6.4f ",
1780 iSt, fMapStationRpcId[iSt], iH, pHit->GetZ(), dDX, dDY, dDZ, dDT, dDTB, dZZ, dTt)
1781 << tPar->ToString();
1782
1783 vhPullX[iSt]->Fill(dDX);
1784 vhPullY[iSt]->Fill(dDY);
1785 vhPullZ[iSt]->Fill(dZZ);
1786 vhPullT[iSt]->Fill(dDT);
1787 vhPullTB[iSt]->Fill(dDTB);
1788 vhTrefRms[iSt]->Fill(dTexpErr);
1789 vhFitDT0[iSt]->Fill(pHit->GetTime() - pTrk->GetT0());
1790 vhFitT0Err[iSt]->Fill(pTrk->GetT0Err());
1791 vhFitTt[iSt]->Fill(pTrk->GetTt());
1792 vhFitTtErr[iSt]->Fill(pTrk->GetTtErr() / pTrk->GetTt());
1793 vhFitDTMean[iSt]->Fill(fTrackletTools->GetDTMean(pTrk, pHit));
1794 vhFitDTMeanErr[iSt]->Fill(fTrackletTools->GetDTMeanError(pTrk, pHit));
1795 vhResidualTBWalk[iSt]->Fill(dTOT, dDTB);
1796 vhResidualYWalk[iSt]->Fill(dTOT, dDY);
1797
1798 fhPullT_Smt->Fill((Double_t) fMapRpcIdParInd[fMapStationRpcId[iSt]], dDT);
1799 fhPullX_Smt->Fill((Double_t) fMapRpcIdParInd[fMapStationRpcId[iSt]], dDX);
1800 fhPullY_Smt->Fill((Double_t) fMapRpcIdParInd[fMapStationRpcId[iSt]], dDY);
1801 bool bCalModule = kTRUE;
1802 if (bCalModule) { /* fill all counters of module*/
1803 int iModType = CbmTofAddress::GetSmType(fMapStationRpcId[iSt]);
1804 if (iModType < 3) {
1807 for (int i = 0; i < fDigiBdfPar->GetNbDet(); i++) {
1808 if (i != iRpc) {
1809 int iRpcId = CbmTofAddress::GetUniqueAddress(iSm, iRpc, 0, 0, iModType);
1810 fhPullX_Smt->Fill((Double_t) fMapRpcIdParInd[iRpcId], dDX);
1811 fhPullY_Smt->Fill((Double_t) fMapRpcIdParInd[iRpcId], dDY);
1812 }
1813 }
1814 }
1815 }
1816 /*
1817 fhPullT_Smt->Fill((Double_t)fMapRpcIdParInd[fMapStationRpcId[iSt]], fTrackletTools->GetTdif(pTrk,fMapStationRpcId[iSt], pHit) );
1818 fhPullX_Smt->Fill((Double_t)fMapRpcIdParInd[fMapStationRpcId[iSt]], fTrackletTools->GetXdif(pTrk,fMapStationRpcId[iSt], pHit) );
1819 fhPullY_Smt->Fill((Double_t)fMapRpcIdParInd[fMapStationRpcId[iSt]], fTrackletTools->GetYdif(pTrk,fMapStationRpcId[iSt], pHit) );
1820 */
1821 fhPullZ_Smt->Fill((Double_t) fMapRpcIdParInd[fMapStationRpcId[iSt]], dZZ);
1822
1823 Double_t dDeltaTt = dTt - fTtTarg;
1824 fhDeltaTt_Smt->Fill((Double_t) fMapRpcIdParInd[fMapStationRpcId[iSt]], dDeltaTt);
1825 if (fiBeamCounter != -1) {
1826 Double_t dDeltaTc = pHit->GetTime() - pRefHit->GetTime() - pHit->GetR() / 29.979;
1827 fhDeltaTc_Smt->Fill((Double_t) fMapRpcIdParInd[fMapStationRpcId[iSt]], dDeltaTc);
1828 }
1829 //XXX use BRef as Referenz!!!
1830 if (pRefHit != NULL) {
1831 if (pHit != pRefHit) {
1832 Double_t dTOff = dDeltaTt * //pHit->GetR();
1833 TMath::Sqrt(TMath::Power(pHit->GetX() - pRefHit->GetX(), 2)
1834 + TMath::Power(pHit->GetY() - pRefHit->GetY(), 2)
1835 + TMath::Power(pHit->GetZ() - pRefHit->GetZ(), 2))
1836 * TMath::Sign(1, pHit->GetZ() - pRefHit->GetZ());
1837 fhTOff_Smt->Fill((Double_t) fMapRpcIdParInd[fMapStationRpcId[iSt]], dTOff);
1838 }
1839 }
1840 } // station loop end
1841
1842 } // condition on NofHits>2 end
1843 Double_t dTt = fTrackletTools->FitTt(pTrk, -1); // restore full fit
1844 if (dTt == 0) {
1845 LOG(error) << "Invalid inverse velocity ";
1846 continue;
1847 }
1848 // monitoring of tracklet hits with selected velocities from reference counters
1849 if (TMath::Abs(pTrk->GetRefVel((UInt_t)(fNReqStations - 1)) - fdRefVelMean) < fdRefDVel) {
1850 fhTrklVelHMul->Fill(pTrk->GetNofHits(), 1. / pTrk->GetTt());
1851 for (Int_t iH = 0; iH < pTrk->GetNofHits(); iH++) {
1852 CbmTofHit* pHit = pTrk->GetTofHitPointer(iH);
1853 Int_t iChId = pHit->GetAddress();
1854 CbmTofCell* fChannelInfo = fDigiPar->GetCell(iChId);
1855 if (NULL == fChannelInfo) continue;
1856 Int_t iAddr = iChId & DetMask;
1857 Int_t iSt = GetStationOfAddr(iAddr);
1858 Double_t hitpos[3] = {3 * 0.};
1859 Double_t hitpos_local[3] = {3 * 0.};
1860 gGeoManager->FindNode(fChannelInfo->GetX(), fChannelInfo->GetY(), fChannelInfo->GetZ());
1861 Int_t iRpcInd = fMapRpcIdParInd[iChId & DetMask];
1862 hitpos[0] = pHit->GetX() - (Double_t) fhPullX_Smt_Off->GetBinContent(iRpcInd + 1);
1863 hitpos[1] = pHit->GetY() - (Double_t) fhPullY_Smt_Off->GetBinContent(iRpcInd + 1);
1864 hitpos[2] = pHit->GetZ() - (Double_t) fhPullZ_Smt_Off->GetBinContent(iRpcInd + 1);
1865 gGeoManager->MasterToLocal(hitpos, hitpos_local);
1866 vhXY_AllTracks[iSt]->Fill(hitpos_local[0], hitpos_local[1]);
1867 }
1868
1869 if (pTrk->GetNofHits() >= fNReqStations) { // all possible hits are there
1870 LOG(debug) << "Complete Tracklet in event " << fiEvent;
1871
1872 for (Int_t iSt = 0; iSt < fNTofStations; iSt++) {
1873 Int_t iH = pTrk->GetStationHitIndex(fMapStationRpcId[iSt]); // Station Hit index
1874 if (iH < 0) {
1875 LOG(debug1) << " Incomplete Tracklet, skip station " << iSt;
1876 continue; // Station not part of tracklet
1877 }
1878 CbmTofHit* pHit = (CbmTofHit*) fTofHitArray->At(iH);
1879 Int_t iChId = pHit->GetAddress();
1880 Double_t hitpos[3] = {3 * 0.};
1881 Double_t hitpos_local[3] = {3 * 0.};
1882 CbmTofCell* fChannelInfo = fDigiPar->GetCell(iChId);
1883 if (NULL == fChannelInfo) {
1884 //faked hit, take init values
1885 }
1886 else {
1887 /* TGeoNode *fNode=*/ // prepare global->local trafo
1888 gGeoManager->FindNode(fChannelInfo->GetX(), fChannelInfo->GetY(), fChannelInfo->GetZ());
1889 Int_t iRpcInd = fMapRpcIdParInd[iChId & DetMask];
1890 hitpos[0] = pHit->GetX() - (Double_t) fhPullX_Smt_Off->GetBinContent(iRpcInd + 1);
1891 hitpos[1] = pHit->GetY() - (Double_t) fhPullY_Smt_Off->GetBinContent(iRpcInd + 1);
1892 hitpos[2] = pHit->GetZ() - (Double_t) fhPullZ_Smt_Off->GetBinContent(iRpcInd + 1);
1893 /* TGeoNode* cNode= */ gGeoManager->GetCurrentNode();
1894 gGeoManager->MasterToLocal(hitpos, hitpos_local);
1895 }
1896 vhXY_AllStations[iSt]->Fill(hitpos_local[0], hitpos_local[1]);
1897 Double_t dDX = pHit->GetX() - pTrk->GetFitX(pHit->GetZ()); // - tPar->GetX() - tPar->GetTx()*dDZ;
1898 Double_t dDY = pHit->GetY() - pTrk->GetFitY(pHit->GetZ()); // - tPar->GetTy()*dDZ;
1899 //Double_t dDT = pHit->GetTime() - pTrk->GetFitT(pHit->GetR()); //pTrk->GetTdif(fMapStationRpcId[iSt]);
1900 Double_t dDTB = fTrackletTools->GetTdif(pTrk, fMapStationRpcId[iSt],
1901 pHit); // ignore pHit in calc of reference
1902 vhXY_DX[iSt]->Fill(hitpos_local[0], hitpos_local[1], dDX);
1903 vhXY_DY[iSt]->Fill(hitpos_local[0], hitpos_local[1], dDY);
1904 vhXY_DT[iSt]->Fill(hitpos_local[0], hitpos_local[1], dDTB);
1905 Double_t dCSZ = ((Double_t)(pHit->GetFlag() % 100)) * 0.5;
1906 Double_t dTOT = ((Double_t) pHit->GetCh()) * 0.1 / dCSZ; // counteract UHIT flagging
1907 vhXY_TOT[iSt]->Fill(hitpos_local[0], hitpos_local[1], dTOT);
1908 vhXY_CSZ[iSt]->Fill(hitpos_local[0], hitpos_local[1], dCSZ);
1909
1910 // Refit tracklet without Dut info
1911 fTrackletTools->Line3DFit(pTrk, iChId & DetMask);
1912 hitpos[0] = pTrk->GetFitX(pHit->GetZ());
1913 hitpos[1] = pTrk->GetFitY(pHit->GetZ());
1914 hitpos[2] = pHit->GetZ();
1915 gGeoManager->MasterToLocal(hitpos, hitpos_local);
1916 vhXY_AllFitStations[iSt]->Fill(hitpos_local[0], hitpos_local[1]);
1917 fTrackletTools->Line3DFit(pTrk, -1);
1918
1919 // debugging consistency of geometry transformations ....
1920 if (fair::Logger::Logging(fair::Severity::debug3)) {
1921 if (iSt == fNReqStations - 1) { // treat as if not found
1922 Int_t iAddr = fMapStationRpcId[iSt];
1923 CbmTofCell* fChannelInfoD = fDigiPar->GetCell(iAddr);
1924 Double_t zPos = 0;
1925 Double_t zPosMiss = -1;
1926 Double_t hitposD[3];
1927 Double_t hitpos_localD[3];
1928 Int_t NIter = 5;
1929 Int_t iRpcInd = fMapRpcIdParInd[iAddr];
1931 while (zPos != zPosMiss && 0 < NIter--) {
1932 fChannelInfoD = fDigiPar->GetCell(iAddr);
1933 gGeoManager->FindNode(fChannelInfoD->GetX(), fChannelInfoD->GetY(), fChannelInfoD->GetZ());
1934 zPos = fChannelInfoD->GetZ() + (Double_t) fhPullZ_Smt_Off->GetBinContent(iRpcInd + 1);
1935 hitposD[0] = pTrk->GetFitX(zPos) - (Double_t) fhPullX_Smt_Off->GetBinContent(iRpcInd + 1);
1936 hitposD[1] = pTrk->GetFitY(zPos) - (Double_t) fhPullY_Smt_Off->GetBinContent(iRpcInd + 1);
1937 hitposD[2] = fChannelInfoD->GetZ();
1938 /* TGeoNode* cNode=*/gGeoManager->GetCurrentNode();
1939 gGeoManager->MasterToLocal(hitposD, hitpos_localD);
1940 // Check for consistency of geometry
1941 Int_t iChTrafo = CbmTofAddress::GetChannelId(iAddr);
1942 Int_t iChMiss = hitpos_localD[0] / fChannelInfoD->GetSizex() + (iNbCh - 1) / 2;
1943 if (iChMiss < 0) iChMiss = 0;
1944 if (iChMiss > iNbCh - 1) iChMiss = iNbCh - 1;
1945 assert(fDigiBdfPar);
1946 if (iChMiss > -1 && iChMiss < iNbCh) {
1947 Int_t iAddrMiss =
1949 iChMiss, 0, CbmTofAddress::GetSmType(iAddr));
1950 CbmTofCell* fChannelInfoMiss = fDigiPar->GetCell(iAddrMiss);
1951 if (NULL == fChannelInfoMiss) {
1952 LOG(fatal) << Form("Geo consistency check 0x%08x, 0x%08x failed at "
1953 "St%d, z=%7.2f,%7.2f: iChTrafo %d, Miss %d , "
1954 "xloc %6.2f, dx %4.2f",
1955 iAddr, iAddrMiss, iSt, zPos, zPosMiss, iChTrafo, iChMiss, hitpos_local[0],
1956 fChannelInfo->GetSizex());
1957 }
1958 zPosMiss = fChannelInfoMiss->GetZ() + (Double_t) fhPullZ_Smt_Off->GetBinContent(iRpcInd + 1);
1959 LOG(debug) << Form("Geo consistency check 0x%08x at St%d, z=%7.2f,%7.2f: "
1960 "iChTrafo %d, Miss %d , xloc %6.2f, dx %4.2f",
1961 iAddr, iSt, zPos, zPosMiss, iChTrafo, iChMiss, hitpos_localD[0],
1962 fChannelInfoD->GetSizex());
1963 fChannelInfo = fChannelInfoMiss;
1964 iAddr = iAddrMiss;
1965 }
1966 else
1967 zPosMiss = zPos;
1968 LOG(debug) << Form("Predicted hit in Last Station 0x%08x at local pos x %6.2f, "
1969 "y %6.2f, z %6.2f, cell %p",
1970 iAddr, hitpos_localD[0], hitpos_localD[1], zPos, fChannelInfoD);
1971 LOG(debug) << Form("Measured hit in Last Station 0x%08x at local pos x %6.2f, "
1972 "y %6.2f, z %6.2f, cell %p",
1973 pHit->GetAddress(), hitpos_local[0], hitpos_local[1], pHit->GetZ(), fChannelInfo);
1974 }
1975 }
1976 }
1977 }
1978 }
1979 else {
1980 if (pTrk->GetNofHits() == fNReqStations - 1) { // one hit missing
1981 for (Int_t iSt = 0; iSt < fNTofStations; iSt++) {
1982 Int_t iH = pTrk->GetStationHitIndex(fMapStationRpcId[iSt]); // Station Hit index
1983 if (iH < 0) { // find geo element for the missing Station iSt
1984 Int_t iAddr = fMapStationRpcId[iSt];
1985 if (iAddr < 1) continue;
1986 //Int_t iChId = CbmTofAddress::GetUniqueAddress(0,0,0,0,iSmType);
1987 //CbmTofCell* fChannelInfo = fDigiPar->GetCell( iChId );
1988 CbmTofCell* fChannelInfo = fDigiPar->GetCell(iAddr);
1989 if (NULL == fChannelInfo) {
1990 if (iWarnNotDefined++ < 100)
1991 LOG(warning) << Form("CbmTofFindTracks::FillHistograms: Cell "
1992 "0x%08x not defined for Station %d",
1993 iAddr, iSt);
1994 continue;
1995 }
1996 /* TGeoNode *fNode= */ // prepare global->local trafo
1997 Double_t zPos = 0;
1998 Double_t zPosMiss = -1;
1999 Double_t hitpos[3];
2000 Double_t hitpos_local[3];
2001 Int_t NIter = 5;
2002 Int_t iRpcInd = fMapRpcIdParInd[iAddr];
2004 while (zPos != zPosMiss && 0 < NIter--) {
2005 fChannelInfo = fDigiPar->GetCell(iAddr);
2006 gGeoManager->FindNode(fChannelInfo->GetX(), fChannelInfo->GetY(), fChannelInfo->GetZ());
2007 zPos = fChannelInfo->GetZ() + (Double_t) fhPullZ_Smt_Off->GetBinContent(iRpcInd + 1);
2008 hitpos[0] = pTrk->GetFitX(zPos) - (Double_t) fhPullX_Smt_Off->GetBinContent(iRpcInd + 1);
2009 hitpos[1] = pTrk->GetFitY(zPos) - (Double_t) fhPullY_Smt_Off->GetBinContent(iRpcInd + 1);
2010 hitpos[2] = fChannelInfo->GetZ();
2011 /* TGeoNode* cNode=*/gGeoManager->GetCurrentNode();
2012 gGeoManager->MasterToLocal(hitpos, hitpos_local);
2013 // Check for consistency of geometry
2014 Int_t iChTrafo = CbmTofAddress::GetChannelId(iAddr);
2015 Int_t iChMiss = hitpos_local[0] / fChannelInfo->GetSizex() + (iNbCh - 1) / 2;
2016 if (iChMiss < 0) iChMiss = 0;
2017 if (iChMiss > iNbCh - 1) iChMiss = iNbCh - 1;
2018 assert(fDigiBdfPar);
2019 if (iChMiss > -1 && iChMiss < iNbCh) {
2020 Int_t iAddrMiss =
2022 iChMiss, 0, CbmTofAddress::GetSmType(iAddr));
2023 CbmTofCell* fChannelInfoMiss = fDigiPar->GetCell(iAddrMiss);
2024 if (NULL == fChannelInfoMiss) {
2025 LOG(fatal) << Form("Geo consistency check 0x%08x, 0x%08x failed at "
2026 "St%d, z=%7.2f,%7.2f: iChTrafo %d, Miss %d , "
2027 "xloc %6.2f, dx %4.2f",
2028 iAddr, iAddrMiss, iSt, zPos, zPosMiss, iChTrafo, iChMiss, hitpos_local[0],
2029 fChannelInfo->GetSizex());
2030 }
2031 zPosMiss = fChannelInfoMiss->GetZ() + (Double_t) fhPullZ_Smt_Off->GetBinContent(iRpcInd + 1);
2032 /*
2033 LOG(debug) << Form("Geo consistency check 0x%08x at St%d, z=%7.2f,%7.2f: "
2034 "iChTrafo %d, Miss %d , xloc %6.2f, dx %4.2f",
2035 iAddr, iSt, zPos, zPosMiss, iChTrafo, iChMiss, hitpos_local[0],
2036 fChannelInfo->GetSizex());
2037 */
2038 fChannelInfo = fChannelInfoMiss;
2039 iAddr = iAddrMiss;
2040 }
2041 else
2042 zPosMiss = zPos;
2043 }
2044 if (iSt == fNReqStations - 1)
2045 LOG(debug) << Form("Missing hit in Last Station in event %d at "
2046 "local pos x %6.2f, y %6.2f, z %6.2f",
2047 fiEvent, hitpos_local[0], hitpos_local[1], zPos);
2048
2049 vhXY_MissedStation[iSt]->Fill(hitpos_local[0], hitpos_local[1]);
2050
2051 // correlation analysis
2052 for (Int_t iTrk1 = iTrk + 1; iTrk1 < fTrackArray->GetEntriesFast(); iTrk1++) {
2053 CbmTofTracklet* pTrk1 = (CbmTofTracklet*) fTrackArray->At(iTrk1);
2054 if (NULL == pTrk1 || pTrk == pTrk1) continue;
2055 if (pTrk1->GetNofHits() >= fNReqStations) { // all possible hits are there
2056 fhTrklDT0SmMis->Fill(iSt, pTrk->GetFitT(0.) - pTrk1->GetFitT(0.));
2057 }
2058 else {
2059 if (pTrk1->GetNofHits() == fNReqStations - 1) { // one hit missing
2060 for (Int_t iSt1 = 0; iSt1 < fNTofStations; iSt1++) {
2061 Int_t iH1 = pTrk1->GetStationHitIndex(fMapStationRpcId[iSt1]); // Station Hit index
2062 if (iH1 < 0) { // find geo element for the missing Station iSt
2063 //Int_t iSmType1 = fMapStationRpcId[iSt1]; (FU) not used
2064 //if (iSmType1 < 1) continue;
2065 fhTrklDT0StMis2->Fill(Double_t(iSt * 10 + iSt1), pTrk->GetFitT(0.) - pTrk1->GetFitT(0.));
2066 }
2067 }
2068 }
2069 }
2070 }
2071 }
2072 }
2073 }
2074 }
2075 } // velocity selection end
2076 } // loop over tracklets end
2077
2078 if (HMul.size() > 3) {
2079 fhTrklMulMaxMM->Fill(HMul[fNReqStations], HMul[fNReqStations - 1]);
2080 }
2081 if (HMul.size() > 5) fhTrklMul3D->Fill(HMul[fNTofStations], HMul[fNTofStations - 1], HMul[fNTofStations - 2]);
2082 fhTrklMulNhits->Fill(fTofHitArray->GetEntriesFast(), iTMul);
2083 fhTrackingTimeNhits->Fill(fTofHitArray->GetEntriesFast(), fdTrackingTime);
2084
2085 // print info about special events
2086 if (0)
2087 if (5 < fNTofStations)
2088 if (HMul[6] > 1) { // temporary
2089 //if (HMul[fNTofStations]>0)
2090 //LOG(info)<<"Found "<<HMul[fNTofStations]<<" max length tracklets in event "<<fiEvent
2091 LOG(info) << "Found " << HMul[6] << " max length tracklets in event " << fiEvent << " within "
2092 << fTofHitArray->GetEntriesFast() << " hits ";
2093 for (Int_t iTrk = 0; iTrk < fTrackArray->GetEntriesFast(); iTrk++) {
2094 CbmTofTracklet* pTrk = (CbmTofTracklet*) fTrackArray->At(iTrk);
2095 if (NULL == pTrk) continue;
2096 pTrk->PrintInfo();
2097 }
2098 }
2099 if (kFALSE) // print event info for special events
2100 if (fTrackArray->GetEntriesFast() > 25) { // temporary
2101 LOG(info) << "Found high track multiplicity of " << fTrackArray->GetEntriesFast() << " in event " << fiEvent
2102 << " from " << fTofHitArray->GetEntriesFast() << " hits ";
2103 for (Int_t iTrk = 0; iTrk < fTrackArray->GetEntriesFast(); iTrk++) {
2104 CbmTofTracklet* pTrk = (CbmTofTracklet*) fTrackArray->At(iTrk);
2105 if (NULL == pTrk) continue;
2106 pTrk->PrintInfo();
2107 }
2108 }
2109
2110 if (iTMul > 1) {
2111 LOG(debug) << Form("CbmTofFindTracks::FillHistograms NTrkl %d(%d) in event %d", iTMul,
2112 fTrackArray->GetEntriesFast(), fiEvent);
2113 for (Int_t iTrk = 0; iTrk < fTrackArray->GetEntriesFast(); iTrk++) {
2114 CbmTofTracklet* pTrk = (CbmTofTracklet*) fTrackArray->At(iTrk);
2115 if (NULL == pTrk) continue;
2116 if (pTrk->GetNofHits() > fMinNofHits) { // for further analysis request min # of matched hits
2117 for (Int_t iTrk1 = iTrk + 1; iTrk1 < fTrackArray->GetEntriesFast(); iTrk1++) {
2118 CbmTofTracklet* pTrk1 = (CbmTofTracklet*) fTrackArray->At(iTrk1);
2119 if (NULL == pTrk1) continue;
2120 if (pTrk1->GetNofHits() > fMinNofHits) { // for further analysis request min # of matched hits
2121 //LOG(info) << " -D- iT "<<iTrk<<", iT1 "<<iTrk1;
2122 fhTrklT0Mul->Fill(iTMul, pTrk->GetFitT(0.) - pTrk1->GetFitT(0.));
2123 }
2124 }
2125 }
2126 }
2127 }
2128
2129 LOG(debug1) << Form("CbmTofFindTracks::FillHistograms: HMul.size() %u ", (UInt_t) HMul.size());
2130 for (UInt_t uHMul = 2; uHMul < HMul.size(); uHMul++) {
2131 LOG(debug1) << Form("CbmTofFindTracks::FillHistograms() HMul %u, #%d", uHMul, HMul[uHMul]);
2132 if (HMul[uHMul] > 0) {
2133 fhTrklHMul->Fill(uHMul, HMul[uHMul]);
2134 }
2135 }
2136
2137 for (Int_t iHit = 0; iHit < fTofHitArray->GetEntriesFast(); iHit++) { // loop over Hits
2138 CbmTofHit* pHit = (CbmTofHit*) fTofHitArray->At(iHit);
2139 // Int_t iSmType = CbmTofAddress::GetSmType( pHit->GetAddress() & DetMask ); (FU) not used
2140 Int_t iAddr = (pHit->GetAddress() & DetMask);
2141 fhAllHitsSmTypes->Fill(GetStationOfAddr(iAddr));
2142 //LOG(info) << " -D- " << iSmType <<", " << fTypeStation[iSmType];
2143 if (GetStationOfAddr(iAddr) > -1) fhAllHitsStation->Fill(GetStationOfAddr(iAddr));
2144 }
2145 // vertex stuff
2146 fhVTXNorm->Fill(fVTXNorm);
2147 if (fVTXNorm > 0.) {
2148 fhVTX_XY0->Fill(fVTX_X, fVTX_Y);
2149 for (Int_t iTrk = 0; iTrk < fTrackArray->GetEntriesFast(); iTrk++) {
2150 CbmTofTracklet* pTrk = (CbmTofTracklet*) fTrackArray->At(iTrk);
2151 if (NULL == pTrk) continue;
2152 if (Double_t w = pTrk->GetNofHits() > (Double_t) fMinNofHits) {
2153 if (fVTXNorm > w) {
2154 Double_t DeltaT0 = pTrk->GetFitT(0.) - (fVTXNorm * fVTX_T - w * pTrk->GetFitT(0.)) / (fVTXNorm - w);
2155 fhVTX_DT0_Norm->Fill(fVTXNorm, DeltaT0);
2156 }
2157 }
2158 }
2159 }
2160 if (0 == fMapStationRpcId[0]) { // Generated Pseudo TofHit at origin
2161 fTofHitArray->RemoveAt(fTofHitArray->GetEntriesFast() - 1); // remove added hit
2162 }
2163}
2164
2166{
2168 for (Int_t i = 0; i < 10; i++)
2169 fTypeStation[i] = -1; // initialize
2170 for (Int_t i = 0; i < fNTofStations; i++) {
2171 Int_t iSm = ival % 10;
2172 Int_t iSt = fNTofStations - 1 - i;
2173 Int_t iAddr = CbmTofAddress::GetUniqueAddress(0, 0, 0, 0, iSm);
2174 fStationType[iSt] = iSm;
2175 fTypeStation[iSm] = iSt;
2176 fMapStationRpcId[iSt] = iAddr;
2177 ival = (ival - iSm) / 10;
2178 }
2179}
2180
2181void CbmTofFindTracks::SetStation(Int_t iVal, Int_t iModType, Int_t iModId, Int_t iRpcId)
2182{
2183 Int_t iCenterCh = 0;
2184 if (NULL != fDigiBdfPar) iCenterCh = TMath::Floor((fDigiBdfPar->GetNbChan(iModType, iRpcId) - 1) / 2);
2185 Int_t iAddr = CbmTofAddress::GetUniqueAddress(iModId, iRpcId, iCenterCh, 0, iModType);
2186 fMapStationRpcId[iVal] = iAddr;
2187 LOG(info) << "SetStation: added " << iVal
2188 << Form(" TSRC %d%d%d%02d, addr 0x%08x ", iModType, iModId, iRpcId, iCenterCh, iAddr);
2189}
2190
2191void CbmTofFindTracks::SetBeamCounter(Int_t iModType, Int_t iModId, Int_t iRpcId)
2192{
2193 fiBeamCounter = CbmTofAddress::GetUniqueAddress(iModId, iRpcId, 0, 0, iModType);
2194}
2195
2197{
2198 std::map<Int_t, Int_t>::iterator it;
2199 Int_t iNSt = 0;
2200 Int_t iSt = -1;
2201 for (it = fMapStationRpcId.begin(); it != fMapStationRpcId.end(); it++) {
2202 //std::map <Int_t, Int_t>::iterator it = fMapStationRpcId.find(iAddr);
2203 /*
2204 Int_t iDetId=it->second;
2205 LOG(debug)<<Form("GetStationOfAddr 0x%08x, TSR %d%d%d: St %d, TSR %d%d%d ", iAddr,
2206 CbmTofAddress::GetSmType(iAddr),CbmTofAddress::GetSmId(iAddr),CbmTofAddress::GetRpcId(iAddr),
2207 it->first, CbmTofAddress::GetSmType(iDetId),CbmTofAddress::GetSmId(iDetId),CbmTofAddress::GetRpcId(iDetId));
2208 */
2209 if (iNSt > fNTofStations) LOG(fatal) << "Invalid NSt " << iNSt;
2210 iNSt++;
2211 if (it->second == iAddr) break;
2212 }
2213 if (it != fMapStationRpcId.end()) iSt = it->first;
2214 /*
2215 if(it->first == fMapStationRpcId.size())
2216 {
2217 PrintSetup();
2218 LOG(fatal)<<Form("CbmTofFindTracks::GetStationOfAddr failed for 0x%08x, found Station = %d",iAddr,it->first)
2219 ;
2220 }
2221 */
2222 return iSt;
2223}
2224
2226{
2227 for (std::map<Int_t, Int_t>::iterator it = fMapStationRpcId.begin(); it != fMapStationRpcId.end(); it++) {
2228 Int_t iDetId = it->second;
2229 LOG(debug) << " <I> Tracking station " << it->first << " contains RpcId "
2230 << Form("0x%08x, TSR %d%d%d", it->second, CbmTofAddress::GetSmType(iDetId),
2232 }
2233 if ((Int_t) fMapStationRpcId.size() > fNTofStations)
2234 LOG(fatal) << "Invalid NTofStations " << fNTofStations << ", " << fMapStationRpcId.size();
2235}
2236
2237Double_t CbmTofFindTracks::GetTOff(Int_t iAddr)
2238{
2239 //LOG(info) << Form(" <D> GetTOff for 0x%08x at HistoIndex %d: %7.1f ", iAddr, fMapRpcIdParInd[iAddr],
2240 //(Double_t)fhPullT_Smt_Off->GetBinContent( fMapRpcIdParInd[iAddr] + 1));
2241 return (Double_t) fhPullT_Smt_Off->GetBinContent(fMapRpcIdParInd[iAddr] + 1);
2242}
2243
2244Double_t CbmTofFindTracks::GetSigT(Int_t iAddr) { return fvTsig[GetStationOfAddr(iAddr)]; }
2245
2246Double_t CbmTofFindTracks::GetSigX(Int_t iAddr) { return fvXsig[GetStationOfAddr(iAddr)]; }
2247
2248Double_t CbmTofFindTracks::GetSigY(Int_t iAddr) { return fvYsig[GetStationOfAddr(iAddr)]; }
2249
2250Double_t CbmTofFindTracks::GetSigZ(Int_t iAddr) { return fvZsig[GetStationOfAddr(iAddr)]; }
2251
2252Double_t CbmTofFindTracks::GetStationSigT(Int_t iSt) { return fvTsig[iSt]; }
2253
2254Double_t CbmTofFindTracks::GetStationSigX(Int_t iSt) { return fvXsig[iSt]; }
2255
2256Double_t CbmTofFindTracks::GetStationSigY(Int_t iSt) { return fvYsig[iSt]; }
2257
2258Double_t CbmTofFindTracks::GetStationSigZ(Int_t iSt) { return fvZsig[iSt]; }
2259
2261{
2262 Int_t iNSt = 0;
2263 for (Int_t iSt = 0; iSt < fNTofStations; iSt++) {
2264 if (fStationHMul[iSt] > 0 && fStationHMul[iSt] < fiStationMaxHMul) iNSt++;
2265 LOG(debug1) << Form("Station %d, Addr 0x%08x, HMul %d, Max %d, Sum %d", iSt, GetAddrOfStation(iSt),
2266 fStationHMul[iSt], fiStationMaxHMul, iNSt);
2267 }
2268 return iNSt;
2269}
2270
2272{
2273 for (Int_t iSt = 0; iSt < fNTofStations; iSt++)
2274 fStationHMul[iSt] = 0;
2275}
2276
2278{
2279 // collect unused hits in active tracking stations
2280 Int_t iNbUHits = 0;
2281 for (Int_t iHit = 0; iHit < fTofHitArray->GetEntriesFast(); iHit++) {
2282 CbmTofHit* pHit = (CbmTofHit*) fTofHitArray->At(iHit);
2283 Int_t iAddr = (pHit->GetAddress() & DetMask);
2284 if (pHit->GetFlag() < 100. && GetStationOfAddr(iAddr) < fNTofStations) {
2285 if (!CheckHit2Track(pHit)) // check whether hit could belong to any track
2286 new ((*fTofUHitArray)[iNbUHits++]) CbmTofHit(*pHit);
2287 }
2288 }
2289}
2290
2292{
2293 Int_t iAddr = (pHit->GetAddress() & DetMask);
2294 Int_t iSt = GetStationOfAddr(iAddr);
2295 if (iSt < 0 || iSt >= GetNofStations()) return kFALSE;
2296 for (Int_t iTrk = 0; iTrk < fTrackArray->GetEntriesFast(); iTrk++) {
2297 CbmTofTracklet* pTrk = (CbmTofTracklet*) fTrackArray->At(iTrk);
2298 if (NULL == pTrk) continue;
2299 Double_t dDX = pHit->GetX() - pTrk->GetFitX(pHit->GetZ());
2300 Double_t dDY = pHit->GetY() - pTrk->GetFitY(pHit->GetZ());
2301 Double_t dDT = pHit->GetTime() - pTrk->GetFitT(pHit->GetZ());
2302 LOG(debug) << Form("Test Trk %d with HMul %d for Addr 0x%08x in station %d "
2303 "with dx %5.1f, dy %5.1f, dt %5.1f",
2304 iTrk, pTrk->GetNofHits(), iAddr, iSt, dDX, dDY, dDT);
2305 if (!(pTrk->ContainsAddr(iAddr))) {
2306 LOG(debug3) << "Fill histo " << vhUDXDY_DT[iSt]->GetName();
2307 vhUDXDY_DT[iSt]->Fill(dDX, dDY, dDT);
2308 }
2309 else {
2310 vhUCDXDY_DT[iSt]->Fill(dDX, dDY, dDT);
2311 }
2312 }
2313 return kFALSE;
2314}
2315
2317{
2318 //fInspectEvent = kTRUE;
2319 for (Int_t iSt = 0; iSt < fNTofStations; iSt++) {
2320 if (fStationHMul[iSt] > fiStationMaxHMul) {
2321 fInspectEvent = kFALSE;
2322 }
2323 else {
2325 fInspectEvent = kFALSE;
2326 }
2327 }
2328 }
2329}
2330
2332{
2333 std::map<Int_t, Int_t>::iterator it = fMapRpcIdParInd.begin();
2334 while (it != fMapRpcIdParInd.end()) {
2335 LOG(info) << Form("MapRpcIdParInd: %d, 0x%08x ", it->second, it->first);
2336 it++;
2337 }
2338}
2339
2341
2342bool CbmTofFindTracks::EvalDoublets(int iI0, int iI1, int iI2, double* dTshift)
2343{
2344 //LOG(info)<<"Evaluate time offsets from 3 doublets in triplet "<<iI0<<" "<<iI1<<" "<<iI2;
2345 //for (int i=0; i<3; i++) dTshift[i]=(double)i;
2346 //LOG(info)<<"Initial Tshift " <<dTshift[0]<<" "<<dTshift[1]<<" "<<dTshift[2];
2347 //return kTRUE;
2348
2349 const double c = 30.; //speed of light in cm/ns
2350 const double dNTrks = 100.;
2351 double A[3] = {3 * 0.};
2352 double D[3] = {3 * 0.};
2353 TString hnameDT[3] = {Form("hDoubletDt_%02d%02d", iI0, iI1), Form("hDoubletDt_%02d%02d", iI0, iI2),
2354 Form("hDoubletDt_%02d%02d", iI1, iI2)};
2355 TString hnameDD[3] = {Form("hDoubletDd_%02d%02d", iI0, iI1), Form("hDoubletDd_%02d%02d", iI0, iI2),
2356 Form("hDoubletDd_%02d%02d", iI1, iI2)};
2357 TString hnameV[3] = {Form("hDoubletV_%02d%02d", iI0, iI1), Form("hDoubletV_%02d%02d", iI0, iI2),
2358 Form("hDoubletV_%02d%02d", iI1, iI2)};
2359 TH1* hDT[3];
2360 TH1* hDD[3];
2361 TH1* hV[3];
2362 for (int i = 0; i < 3; i++) {
2363 hDT[i] = (TH1*) gROOT->FindObjectAny(hnameDT[i]);
2364 if (NULL == hDT[i]) {
2365 LOG(warn) << "CalHisto " << hnameDT[i] << " not existing";
2366 return kFALSE;
2367 }
2368 hDD[i] = (TH1*) gROOT->FindObjectAny(hnameDD[i]);
2369 if (NULL == hDD[i]) {
2370 LOG(warn) << "CalHisto " << hnameDD[i] << " not existing";
2371 return kFALSE;
2372 }
2373 hV[i] = (TH1*) gROOT->FindObjectAny(hnameV[i]);
2374 if (NULL == hV[i]) {
2375 LOG(warn) << "CalHisto " << hnameV[i] << " not existing";
2376 return kFALSE;
2377 }
2378 if (hDT[i]->Integral() < dNTrks) {
2379 LOG(warn) << "Too few entries for triplet " << iI0 << " " << iI1 << " " << iI2 << " in " << hnameV[i] << ": "
2380 << hV[i]->Integral();
2381 return kFALSE; // return without shifts
2382 }
2383 //D[i]=hDT[i]->GetMean()*hV[i]->GetMean();
2384 D[i] = hDD[i]->GetMean();
2385 A[i] = (D[i] - c * hDT[i]->GetMean()) / c;
2386 }
2387 dTshift[0] = 0.5 * (A[0] + A[1] - A[2]);
2388 dTshift[1] = 0.5 * (-A[0] + A[1] - A[2]);
2389 dTshift[2] = 0.5 * (A[0] + A[1] + A[2]);
2390 LOG(info) << "Final Tshift for " << iI0 << " " << iI1 << " " << iI2 << ": " << dTshift[0] << " " << dTshift[1] << " "
2391 << dTshift[2];
2392 return kTRUE;
2393}
2394
2396{
2397 if ((int) fiStationStatus.size() < iStation) {
2398 fiStationStatus.resize(iStation + 1);
2399 for (int i = 0; i < iStation + 1; i++)
2400 fiStationStatus[i] = 0;
2401 LOG(info) << "StationStatus initialized for " << fiStationStatus.size() << " stations";
2402 }
2403 fiStationStatus[iStation] = iStatus;
2404}
2406{
2407 if (iStation < (int) fiStationStatus.size())
2408 return fiStationStatus[iStation];
2409 else
2410 return 0;
2411}
2412
2413int CbmTofFindTracks::GetNbHits() { return fTofHitArray->GetEntriesFast(); } // for internal use by Calibrator
2414
2416
2418{
2419 if ((int) fTofHitIndexArray.size() > iHit) {
2420 //LOG(debug) << "Hit " << iHit << ", index " << fTofHitIndexArray[iHit];
2421 return fTofHitIndexArray[iHit];
2422 }
2423 else
2424 return iHit;
2425}
static TFile * fHist
TClonesArray * fEventsColl
const Int_t DetMask
const Int_t DetMask
ClassImp(CbmTofFindTracks)
static int iTS
static Int_t iWarnNotDefined
@ k12b
@ k21a
@ k14a
Class characterising one event by a collection of links (indices) to data objects,...
Definition CbmEvent.h:34
size_t GetNofData() const
Definition CbmEvent.cxx:53
uint32_t GetIndex(ECbmDataType type, uint32_t iData)
Definition CbmEvent.cxx:42
void AddData(ECbmDataType type, uint32_t index)
Definition CbmEvent.cxx:33
void SetTimeError(double error)
Definition CbmHit.h:91
double GetTimeError() const
Definition CbmHit.h:77
double GetDz() const
Definition CbmHit.h:72
double GetTime() const
Definition CbmHit.h:76
int32_t GetAddress() const
Definition CbmHit.h:74
void SetZ(double z)
Definition CbmHit.h:80
double GetZ() const
Definition CbmHit.h:71
void SetTime(double time)
Definition CbmHit.h:85
double GetDy() const
Definition CbmPixelHit.h:76
double GetDx() const
Definition CbmPixelHit.h:75
void SetX(double x)
Definition CbmPixelHit.h:92
double GetY() const
Definition CbmPixelHit.h:74
void SetY(double y)
Definition CbmPixelHit.h:93
void SetPositionError(const TVector3 &dpos)
double GetX() const
Definition CbmPixelHit.h:73
void SetPosition(const TVector3 &pos)
Sets position of the hit.
static uint32_t GetUniqueAddress(uint32_t Sm, uint32_t Rpc, uint32_t Channel, uint32_t Side=0, uint32_t SmType=0, uint32_t RpcType=0)
static int32_t GetSmId(uint32_t address)
static int32_t GetRpcId(uint32_t address)
static int32_t GetSmType(uint32_t address)
static int32_t GetChannelId(uint32_t address)
contains filling and updating of calibration histos
void FillCalHist(CbmTofHit *pHit, Int_t iOpt, CbmEvent *pEvent=NULL)
void SetR0Lim(Double_t dVal)
void SetBeam(Bool_t bVal)
Bool_t UpdateCalHist(Int_t iOpt)
Double_t GetSizey() const
Definition CbmTofCell.h:40
Double_t GetY() const
Definition CbmTofCell.h:36
Double_t GetSizex() const
Definition CbmTofCell.h:39
Double_t GetX() const
Definition CbmTofCell.h:35
Double_t GetZ() const
Definition CbmTofCell.h:37
virtual InitStatus Init()
virtual int32_t GetSModule(const int32_t detectorId)=0
virtual int32_t GetCounter(const int32_t detectorId)=0
virtual int32_t GetCell(const int32_t detectorId)=0
virtual int32_t GetSMType(const int32_t detectorId)=0
Parameters class for the CBM ToF digitizer using beam data distributions.
Int_t GetNbChan(Int_t iSmType, Int_t iRpc) const
Int_t GetNbTrackingStations() const
Int_t GetDetUId(Int_t iDet)
Int_t GetTrackingStation(Int_t iSmType, Int_t iSm, Int_t iRpc) const
Double_t GetSigVel(Int_t iSmType, Int_t iSm, Int_t iRpc) const
CbmTofCell * GetCell(Int_t i)
Int_t GetNrOfModules()
Int_t GetCellId(Int_t i)
Int_t fTypeStation[1000]
TClonesArray * fTofUHitArray
std::vector< TH3 * > vhXY_DX
TClonesArray * fTofUHitArrayOut
std::vector< TH3 * > vhXY_TOT
std::vector< Int_t > fRpcAddr
std::vector< TH2 * > vhResidualTBWalk
std::vector< TH3 * > vhXY_CSZ
std::map< Int_t, Int_t > fMapStationRpcId
Double_t GetStationSigZ(Int_t iSt)
std::vector< TH3 * > vhUDXDY_DT
Double_t GetSigZ() const
TClonesArray * fTofHitArrayIn
std::vector< TH1 * > vhPullX
virtual void Exec(Option_t *opt)
CbmTofTrackletTools * fTrackletTools
std::vector< TH1 * > vhTrefRms
std::vector< Double_t > fvTsig
std::vector< Double_t > fvZoff
Double_t GetStationSigT(Int_t iSt)
std::vector< TH2 * > vhXY_AllTracks
Int_t GetStationOfAddr(Int_t iAddr)
std::vector< TH1 * > vhPullTB
TClonesArray * fTofHitArrayOut
std::vector< TH1 * > vhFitTt
TClonesArray * fTrackArrayOut
std::vector< Double_t > fvZsig
Double_t GetSigY() const
std::map< Int_t, Int_t > fMapRpcIdParInd
std::vector< Double_t > fvXsig
std::vector< TH2 * > vhResidualYWalk
void SetBeamCounter(Int_t iModType, Int_t iModId, Int_t iRpcId)
std::vector< Double_t > fvXoff
void SetStations(Int_t ival)
std::vector< Int_t > fStationType
CbmTofDigiPar * fDigiPar
Double_t GetStationSigX(Int_t iSt)
int GetStationStatus(int iStation)
std::vector< TH2 * > vhXY_AllStations
std::vector< TH1 * > vhFitTtErr
virtual void FillHistograms(CbmEvent *tEvent=NULL)
virtual void CreateHistograms()
virtual void FindVertex()
int GetHitIndex(int iHit)
std::vector< int > fiStationStatus
Int_t GetMinNofHits() const
CbmTofHit * GetHitPointer(int iHit)
void SetStationStatus(int iStation, int iStatus)
Double_t GetSigT() const
std::vector< TH1 * > vhPullZ
std::vector< Double_t > fvYoff
virtual Bool_t CheckHit2Track(CbmTofHit *pHit)
std::vector< TH1 * > vhPullT
virtual InitStatus Init()
virtual void SetParContainers()
std::vector< TH1 * > vhFitT0Err
std::vector< Int_t > fStationHMul
std::vector< TH2 * > vhXY_MissedStation
bool EvalDoublets(int iI0, int iI1, int iI2, double *dTshift)
TClonesArray * fTrackArray
std::vector< TH1 * > vhFitDT0
CbmTofTrackFinder * fFinder
void SetStation(Int_t iVal, Int_t iModType, Int_t iModId, Int_t iRpcId)
CbmTofDigiBdfPar * fDigiBdfPar
Double_t GetStationSigY(Int_t iSt)
std::vector< Int_t > fTofHitIndexArray
Double_t GetTOff(Int_t iAddr)
CbmTofGeoHandler * fGeoHandler
TClonesArray * fTofMatchArrayIn
virtual void ExecFind(Option_t *opt, CbmEvent *tEvent=NULL)
CbmTofCalibrator * fTofCalibrator
Double_t GetSigX() const
Int_t GetAddrOfStation(Int_t iVal)
virtual void FillUHits()
Bool_t fbRemoveSignalPropagationTime
std::vector< TH3 * > vhXY_DT
std::vector< TH2 * > vhXY_AllFitStations
Int_t GetNReqStations() const
CbmTofDetectorId * fTofId
static CbmTofFindTracks * fInstance
std::vector< TH1 * > vhFitDTMean
std::vector< TH1 * > vhFitDTMeanErr
std::vector< TH3 * > vhXY_DY
std::vector< TH1 * > vhPullY
std::vector< Double_t > fvYsig
TClonesArray * fTofHitArray
std::vector< TH3 * > vhUCDXDY_DT
TClonesArray * fEventsColl
std::vector< Double_t > fvToff
void MarkStationFired(Int_t iSt)
Int_t Init(Bool_t isSimulation=kFALSE)
int32_t GetCh() const
Definition CbmTofHit.h:76
double GetR() const
Definition CbmTofHit.h:78
int32_t GetFlag() const
Definition CbmTofHit.h:75
virtual Int_t DoFind(TClonesArray *hits, TClonesArray *tracks)=0
virtual void Init()
double GetZy(double Y) const
std::string ToString() const
Return string representation of class.
contains fits and resolution functions
virtual Double_t GetTexpected(CbmTofTracklet *pTrk, Int_t iDetId, CbmTofHit *pHit, double TtLight=0.)
virtual Double_t GetTexpectedError(CbmTofTracklet *pTrk, Int_t iDetId, CbmTofHit *pHit, Double_t dTexp)
Double_t * Line3DFit(CbmTofTracklet *pTrk, Int_t iDetId)
virtual Double_t GetTdif(CbmTofTracklet *pTrk, Int_t iDetId, CbmTofHit *pHit)
virtual Double_t GetDTMeanError(CbmTofTracklet *pTrk, CbmTofHit *pHit)
virtual Double_t GetDTMean(CbmTofTracklet *pTrk, CbmTofHit *pHit)
Double_t FitTt(CbmTofTracklet *pTrk, Int_t iDetId)
Provides information on attaching a TofHit to a TofTrack.
void SetTofHitIndex(int32_t ind, int32_t i)
double GetFitY(double Z)
int32_t GetTofHitIndex(int32_t ind) const
int32_t GetStationHitIndex(int32_t iSm) const
double GetRefVel(uint32_t N)
virtual bool ContainsAddr(int32_t iAddr)
virtual void PrintInfo()
double GetTt() const
int32_t GetNofHits() const
CbmTofHit * GetTofHitPointer(int32_t ind)
double GetFitT(double R)
int32_t GetTofDetIndex(int32_t ind) const
double GetT0Err() const
double GetTtErr() const
double GetT0() const
double GetChiSq() const
CbmTofTrackletParam * GetTrackParameter()
double GetFitX(double Z)