CbmRoot
CbmTofDigiBdfPar.cxx
Go to the documentation of this file.
1 
5 #include "CbmTofDigiBdfPar.h"
6 
7 #include "CbmTofAddress.h" // for CbmTofAddress
8 
9 #include <FairLogger.h> // for LOG, Logger
10 #include <FairParGenericSet.h> // for FairParGenericSet
11 #include <FairParamList.h> // for FairParamList
12 
13 #include <TArrayD.h> // for TArrayD
14 #include <TArrayI.h> // for TArrayI
15 #include <TDirectory.h> // for TDirectory, gDirectory
16 #include <TFile.h> // for TFile
17 #include <TFitResult.h> // for TFitResult
18 #include <TFitResultPtr.h> // for TFitResultPtr
19 #include <TH1.h> // for TH1
20 #include <TH2.h> // for TH2D
21 #include <TMath.h> // for Power, Sqrt
22 #include <TROOT.h> // for TROOT, gROOT
23 
25 
26  CbmTofDigiBdfPar::CbmTofDigiBdfPar(const char* name,
27  const char* title,
28  const char* context)
29  : FairParGenericSet(name, title, context)
30  , fbUseExpDigi(kTRUE)
31  , fbUseOnlyPrim(kFALSE)
32  , fbOneGapTrack(kTRUE)
33  , fiClusterModel(0)
34  , fdFeeGainSigma(0.0)
35  , fdFeeTotThr(0.0)
36  , fdTimeResElec(0.0)
37  , fdTimeResStart(0.0)
38  , fdDeadtime(5.0)
39  , fdSignalPropSpeed(0.0)
40  , fiNbSmTypes(0)
41  , fiNbSm()
42  , fiNbRpc()
43  , fiNbGaps()
44  , fdGapSize()
45  , fdSigVel()
46  , fiNbCh()
47  , fiChType()
48  , fiChOrientation()
49  , fiDetUId()
50  , fMapDetInd()
51  , fsBeamInputFile("")
52  , fsBeamCalibFile("")
53  , fiClusterRadiusModel(-1)
54  , fiSmTypeInpMapp()
55  , fdEfficiency()
56  , fdGapsEfficiency()
57  , fdTimeResolution()
58  , fh1ClusterSize()
59  , fh1ClusterTot()
60  , fdLandauMpv()
61  , fdLandauSigma()
62  , fbMulUseTrackId(kTRUE)
63  , fdMaxTimeDistClust(0.5)
64  , fdMaxSpaceDistClust(2.5) {
65  detName = "Tof";
66 }
67 
69  LOG(DEBUG4) << "Enter CbmTofDigiBdfPar destructor";
70  clear();
71  LOG(DEBUG4) << "Leave CbmTofDigiBdfPar destructor";
72 }
73 
75  status = kFALSE;
76  resetInputVersions();
77  ClearHistos();
78 }
80  for (UInt_t uSmType = 0; uSmType < fh1ClusterSize.size(); uSmType++)
81  delete fh1ClusterSize[uSmType];
82  for (UInt_t uSmType = 0; uSmType < fh1ClusterTot.size(); uSmType++)
83  delete fh1ClusterTot[uSmType];
84  fh1ClusterSize.clear();
85  fh1ClusterTot.clear();
86 }
87 
88 void CbmTofDigiBdfPar::putParams(FairParamList* l) {
89  if (!l) { return; }
90 
91  l->add("UseExpDigi", (Int_t) fbUseExpDigi);
92  l->add("UseOnlyPrim", (Int_t) fbUseOnlyPrim);
93  l->add("ClusterModel", fiClusterModel);
94  l->add("FeeGainSigma", fdFeeGainSigma);
95  l->add("FeeTotThr", fdFeeTotThr);
96  l->add("TimeResElec", fdTimeResElec);
97  l->add("TimeResStart", fdTimeResStart);
98  l->add("Deadtime", fdDeadtime);
99  l->add("SignalPropSpeed", fdSignalPropSpeed);
100  l->add("NbSmTypes", fiNbSmTypes);
101  l->add("NbSm", fiNbSm);
102  l->add("NbRpc", fiNbRpc);
103  // l->add("NbCh", fiNbCh);
104  // l->add("ChType", fiChType);
105  // l->add("ChOrientation", fiChOrientation);
106  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++)
107  l->add(Form("NbGaps%03d", iSmType), fiNbGaps[iSmType]);
108  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++)
109  l->add(Form("GapSize%03d", iSmType), fdGapSize[iSmType]);
110  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++)
111  l->add(Form("SigVel%03d", iSmType), fdSigVel[iSmType]);
112  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++)
113  l->add(Form("NbCh%03d", iSmType), fiNbCh[iSmType]);
114  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++)
115  l->add(Form("ChType%03d", iSmType), fiChType[iSmType]);
116  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++)
117  l->add(Form("ChOrientation%03d", iSmType), fiChOrientation[iSmType]);
118  l->add("BeamInputFile", fsBeamInputFile);
119  l->add("BeamCalibFile", fsBeamCalibFile);
120  l->add("ClusterRadiusModel", fiClusterRadiusModel);
121  l->add("SmTypeInpMapp", fiSmTypeInpMapp);
122  if (1 == fiClusterRadiusModel) {
123  l->add("RadiusLandauMpv", fdLandauMpv);
124  l->add("RadLandauSigma", fdLandauSigma);
125  } // if( 1 == fiClusterRadiusModel )
126  l->add("MulUseTrkId", (Int_t) fbUseOnlyPrim);
127  l->add("MaxTimeDistClust", fdMaxTimeDistClust);
128  l->add("MaxSpaceDistClust", fdMaxSpaceDistClust);
129 }
130 
131 Bool_t CbmTofDigiBdfPar::getParams(FairParamList* l) {
132  if (!l) { return kFALSE; }
133 
134  LOG(debug2) << "Get the tof digitization parameters.";
135 
136  Int_t iTemp;
137  if (!l->fill("UseExpDigi", &iTemp)) return kFALSE;
138  fbUseExpDigi = (1 == iTemp ? kTRUE : kFALSE);
139 
140  if (!l->fill("UseOnlyPrim", &iTemp)) return kFALSE;
141  fbUseOnlyPrim = (1 == iTemp ? kTRUE : kFALSE);
142 
143  if (!l->fill("UseOneGapPerTrk", &iTemp)) return kFALSE;
144  fbOneGapTrack = (1 == iTemp ? kTRUE : kFALSE);
145 
146  if (!l->fill("ClusterModel", &fiClusterModel)) return kFALSE;
147  if (!l->fill("FeeGainSigma", &fdFeeGainSigma)) return kFALSE;
148  if (!l->fill("FeeTotThr", &fdFeeTotThr)) return kFALSE;
149  if (!l->fill("TimeResElec", &fdTimeResElec)) return kFALSE;
150  if (!l->fill("TimeResStart", &fdTimeResStart)) return kFALSE;
151  if (!l->fill("Deadtime", &fdDeadtime))
152  LOG(debug) << "Use default FEE deadtime of " << fdDeadtime << " ns ";
153  if (!l->fill("SignalPropSpeed", &fdSignalPropSpeed)) return kFALSE;
154  if (!l->fill("NbSmTypes", &fiNbSmTypes)) return kFALSE;
155 
156  LOG(debug2) << "CbmTofDigiBdfPar => There are " << fiNbSmTypes
157  << " types of SM to be initialized.";
158 
159 
160  fiNbSm.Set(fiNbSmTypes);
161  if (!l->fill("NbSm", &fiNbSm)) return kFALSE;
162 
163  fiNbRpc.Set(fiNbSmTypes);
164  if (!l->fill("NbRpc", &fiNbRpc)) return kFALSE;
165 
166  fiNbGaps.resize(fiNbSmTypes);
167  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++) {
168  fiNbGaps[iSmType].Set(fiNbRpc[iSmType]);
169  if (!l->fill(Form("NbGaps%03d", iSmType), &(fiNbGaps[iSmType])))
170  return kFALSE;
171  } // for( Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType ++)
172 
173  fdGapSize.resize(fiNbSmTypes);
174  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++) {
175  fdGapSize[iSmType].Set(fiNbRpc[iSmType]);
176  if (!l->fill(Form("GapSize%03d", iSmType), &(fdGapSize[iSmType])))
177  return kFALSE;
178  } // for( Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType ++)
179 
180  fdSigVel.resize(fiNbSmTypes);
181  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++) {
182  fdSigVel[iSmType].Set(fiNbSm[iSmType] * fiNbRpc[iSmType]);
183  if (!l->fill(Form("SigVel%03d", iSmType), &(fdSigVel[iSmType]))) {
184  LOG(WARNING)
185  << "CbmTofDigiBdfPar::getParams => parameter "
186  << Form("SigVel%03d", iSmType) << " not found in the text file. "
187  << "This is normal for geometries < 14a but may indicate file "
188  "corruption "
189  << " for newer geometries. Values are set to default 0.0 cm/ps.";
190  for (Int_t iRpc = 0; iRpc < fiNbRpc[iSmType]; iRpc++)
191  fdSigVel[iSmType].SetAt(0.0, iRpc);
192  //return kFALSE;
193  } // if ( ! l->fill( Form("SigVel%03d", iSmType), &(fdSigVel[iSmType]) ) )
194 
195  if ((fdSigVel[iSmType]).GetSize() < fiNbSm[iSmType] * fiNbRpc[iSmType]) {
196  if ((fdSigVel[iSmType]).GetSize() == fiNbRpc[iSmType]) {
197  LOG(info) << Form(
198  "SigVel%03d Array size: %d < request %d, copy values to all modules ",
199  iSmType,
200  (Int_t)(fdSigVel[iSmType]).GetSize(),
201  fiNbSm[iSmType] * fiNbRpc[iSmType]);
202 
203  TArrayD temp((fdSigVel[iSmType]).GetSize()); // temporary copy of data
204  for (Int_t i = 0; i < (fdSigVel[iSmType]).GetSize(); i++)
205  temp[i] = fdSigVel[iSmType][i];
206 
207  fdSigVel[iSmType].Set(fiNbSm[iSmType] * fiNbRpc[iSmType]);
208  for (Int_t iSm = 0; iSm < fiNbSm[iSmType];
209  iSm++) { // fill elements for all SMs
210  for (Int_t iRpc = 0; iRpc < fiNbRpc[iSmType]; iRpc++)
211  fdSigVel[iSmType][iSm * fiNbRpc[iSmType] + iRpc] = temp[iRpc];
212  }
213  } // if( (fdSigVel[iSmType]).GetSize() == fiNbRpc[iSmType] )
214  else {
215  LOG(error) << "CbmTofDigiBdfPar::getParams => parameter "
216  << Form("SigVel%03d", iSmType)
217  << " has a not matching number of fields: "
218  << (fdSigVel[iSmType]).GetSize() << " instead of "
219  << fiNbRpc[iSmType] << " or "
220  << fiNbSm[iSmType] * fiNbRpc[iSmType];
221  return kFALSE;
222  } // else of if( (fdSigVel[iSmType]).GetSize() == fiNbRpc[iSmType] )
223  }
224 
225  } // for( Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType ++)
226 
227  fiNbCh.resize(fiNbSmTypes);
228  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++) {
229  fiNbCh[iSmType].Set(fiNbRpc[iSmType]);
230  if (!l->fill(Form("NbCh%03d", iSmType), &(fiNbCh[iSmType]))) return kFALSE;
231  } // for( Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType ++)
232 
233  fiChType.resize(fiNbSmTypes);
234  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++) {
235  fiChType[iSmType].Set(fiNbRpc[iSmType]);
236  if (!l->fill(Form("ChType%03d", iSmType), &(fiChType[iSmType])))
237  return kFALSE;
238  } // for( Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType ++)
239 
241  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++) {
242  fiChOrientation[iSmType].Set(fiNbRpc[iSmType]);
243  if (!l->fill(Form("ChOrientation%03d", iSmType),
244  &(fiChOrientation[iSmType])))
245  return kFALSE;
246  } // for( Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType ++)
247 
248 
249  Int_t nDet = 0;
250  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++) {
251  nDet += GetNbSm(iSmType) * GetNbRpc(iSmType);
252  }
253 
254  fiDetUId.Set(nDet); // dimension unique identifier array
255  fMapDetInd.clear(); // reset
256 
257  Int_t iDet = 0;
258  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++) {
259  for (Int_t iSm = 0; iSm < GetNbSm(iSmType); iSm++) {
260  for (Int_t iRpc = 0; iRpc < GetNbRpc(iSmType); iRpc++) {
261  Int_t iAddr = CbmTofAddress::GetUniqueAddress(iSm, iRpc, 0, 0, iSmType);
262  fMapDetInd[iAddr] = iDet;
263  fiDetUId[iDet++] = iAddr;
264  }
265  }
266  } // for( Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType ++)
267  /*
268  fiChType.Set(fiNbSmTypes);
269  if ( ! l->fill("", &) ) return kFALSE;
270 
271  .Set(fiNbSmTypes);
272  if ( ! l->fill("", &) ) return kFALSE;
273 */
274  // Use Text_t as FairParamList does not seem to support TStrings yet
275 
276 
277  Int_t iMaxSizeFilename = 500;
278  Text_t* sTempText;
279  sTempText = new Text_t[iMaxSizeFilename];
280  if (!l->fill("BeamInputFile", sTempText, iMaxSizeFilename)) return kFALSE;
281  fsBeamInputFile = sTempText;
282  if (l->fill("BeamCalibFile", sTempText, iMaxSizeFilename))
283  fsBeamCalibFile = sTempText;
284 
285  if (!l->fill("ClusterRadiusModel", &fiClusterRadiusModel)) return kFALSE;
286 
288  if (!l->fill("SmTypeInpMapp", &fiSmTypeInpMapp)) return kFALSE;
289 
290  if (1 <= fiClusterRadiusModel) {
292  if (!l->fill("RadiusLandauMpv", &fdLandauMpv)) return kFALSE;
294  if (!l->fill("RadLandauSigma", &fdLandauSigma)) return kFALSE;
295  } // if( 1 <= fiClusterRadiusModel )
296 
297  if (!l->fill("MulUseTrkId", &iTemp)) return kFALSE;
298  fbMulUseTrackId = (1 == iTemp ? kTRUE : kFALSE);
299 
300  if (!l->fill("MaxTimeDistClust", &fdMaxTimeDistClust)) return kFALSE;
301  if (!l->fill("MaxSpaceDistClust", &fdMaxSpaceDistClust)) return kFALSE;
302 
303  // Extract beamtime parameters from file
306  // So that the histos are not created in more than 1 exemplar in case of a double call to this method
307  ClearHistos();
308  fh1ClusterSize.resize(fiNbSmTypes);
309  fh1ClusterTot.resize(fiNbSmTypes);
310 
311  return kTRUE;
312 }
313 
315  TDirectory* oldir = gDirectory;
316 
317  TFile* fBeamtimeInput = new TFile(fsBeamInputFile, "READ");
318  if (kFALSE == fBeamtimeInput->IsOpen()) {
319  LOG(error) << "CbmTofDigiBdfPar => Could not open the beamtime data file.";
320  return kFALSE;
321  } // if( kFALSE == fBeamtimeInput->IsOpen() )
322  TArrayD* pInputEff;
323  TArrayD* pInputRes;
324  fBeamtimeInput->GetObject("Efficiency", pInputEff);
325  fBeamtimeInput->GetObject("TimeResol", pInputRes);
326  if (0 == pInputEff) {
327  LOG(error) << "CbmTofDigiBdfPar => Could not recover the Efficiency array "
328  "from the beamtime data file.";
329  return kFALSE;
330  } // if( 0 == pInputEff)
331  if (0 == pInputRes) {
332  LOG(error) << "CbmTofDigiBdfPar => Could not recover the Time Resolution "
333  "array from the beamtime data file.";
334  gDirectory->cd(oldir->GetPath());
335  fBeamtimeInput->Close();
336  return kFALSE;
337  } // if( 0 == pInputEff)
338  if (0 == pInputEff->GetSize() || 0 == pInputRes->GetSize()
339  || pInputEff->GetSize() != pInputRes->GetSize()) {
340  LOG(error) << "CbmTofDigiBdfPar => Efficiency or Time Resolution array "
341  "from the beamtime data file have wrong size.";
342  gDirectory->cd(oldir->GetPath());
343  fBeamtimeInput->Close();
344  return kFALSE;
345  } // if wrong array size
346 
347  TH1* pH1Temp;
348  gROOT->cd();
350  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++)
351  if (fiSmTypeInpMapp[iSmType] < pInputEff->GetSize()) {
352  fdEfficiency[iSmType] = pInputEff->GetAt(fiSmTypeInpMapp[iSmType]);
353 
354  // For now use a simple gap treatment:
355  // - a single fired gap fires the channel
356  // - all gaps have exactly equivalent efficiency/firing probability
357  // - independent gap firing (no charge accumulation or such)
358  // The probability for a hit not to fire at all is then
359  // (1-p)^NGaps with p = gap efficiency and Ngaps the number of gaps in this RPC
360  // This results in the relation: p = 1 - (1 - P)^(1/Ngaps)
361  // with P = RPC efficiency from beamtime
362  fdGapsEfficiency[iSmType].Set(fiNbRpc[iSmType]);
363  fdTimeResolution[iSmType] = pInputRes->GetAt(fiSmTypeInpMapp[iSmType]);
364 
365  for (Int_t iRpc = 0; iRpc < fiNbRpc[iSmType]; iRpc++) {
366  LOG(debug1) << iSmType << " Eff: " << iRpc << " "
367  << fdEfficiency[iSmType] << " " << GetNbGaps(iSmType, iRpc)
368  << " " << (Double_t) GetNbGaps(iSmType, iRpc) << " TRes "
369  << fdTimeResolution[iSmType];
370  fdGapsEfficiency[iSmType][iRpc] =
371  1
372  - TMath::Power(1.0 - fdEfficiency[iSmType],
373  1.0 / (Double_t) GetNbGaps(iSmType, iRpc));
374  }
375 
376 
377  // Cluster Size histograms
378  fBeamtimeInput->GetObject(
379  Form("h1ClusterSizeType%03d", fiSmTypeInpMapp[iSmType]), pH1Temp);
380  if (0 == pH1Temp) {
381  LOG(error) << "CbmTofDigiBdfPar => Could not recover the Cluster Size "
382  "histogram for Sm Type "
383  << iSmType << ", mapped to input type "
384  << fiSmTypeInpMapp[iSmType];
385  gDirectory->cd(oldir->GetPath());
386  fBeamtimeInput->Close();
387  return kFALSE;
388  } // if( 0 == pH1Temp )
389  fh1ClusterSize[iSmType] =
390  (TH1*) (pH1Temp->Clone(Form("ClusterSizeType%03d", iSmType)));
391 
392  // Cluster Total ToT histograms
393  fBeamtimeInput->GetObject(
394  Form("h1ClusterTot%03d", fiSmTypeInpMapp[iSmType]), pH1Temp);
395  if (0 == pH1Temp) {
396  LOG(error) << "CbmTofDigiBdfPar => Could not recover the Cluster Size "
397  "histogram for Sm Type "
398  << iSmType << ", mapped to input type "
399  << fiSmTypeInpMapp[iSmType];
400  gDirectory->cd(oldir->GetPath());
401  fBeamtimeInput->Close();
402  return kFALSE;
403  } // if( 0 == pH1Temp )
404  fh1ClusterTot[iSmType] =
405  (TH1*) (pH1Temp->Clone(Form("ClusterTot%03d", iSmType)));
406  } // if( fiSmTypeInpMapp[iSmType] < pInputEff->GetSize() )
407  else {
408  LOG(error) << "CbmTofDigiBdfPar => Wrong mapping index for Sm Type "
409  << iSmType << ": Out of input boundaries";
410  gDirectory->cd(oldir->GetPath());
411  fBeamtimeInput->Close();
412  return kFALSE;
413  }
414 
415  if (2 == fiClusterRadiusModel) {
417  } // if( 2 == fiClusterRadiusModel )
418  gDirectory->cd(oldir->GetPath());
419  fBeamtimeInput->Close();
420 
421  return kTRUE;
422 }
423 /************************************************************************************/
425  TDirectory* oldir = gDirectory;
426 
427  TFile* fSimInput =
428  new TFile("RadToClustDist_0000_1000_0010_00025_05025_00025.root", "READ");
429  if (kFALSE == fSimInput->IsOpen()) {
430  LOG(error) << "CbmTofDigiBdfPar => Could not open the Simulated Landau "
431  "distribution data file.";
432  return kFALSE;
433  } // if( kFALSE == fSimInput->IsOpen() )
434 
435  // Get the input histograms from the Simulated Landau distributions file
436  // R0 value of the Landau distribution used to generate a simple cluster size distribution as function of
437  // the MPV and sigma values of a Landau fit on said Cluster size distribution
438  TH2D* hFitR0All;
439  fSimInput->GetObject("ClustSizeFitR0All", hFitR0All);
440  // Input Sigma value of the Landau distribution used to generate a simple cluster size distribution as function of
441  // the MPV and sigma values of a Landau fit on said Cluster size distribution
442  TH2D* hFitSigInAll;
443  fSimInput->GetObject("ClustSizeFitSigInAll", hFitSigInAll);
444  // Number of entries for each bin of the R0 of a Landau Input distribution as function of
445  // the MPV and sigma values of a Landau fit on the corresponding Cluster size distribution
446  // The (R0, Input Sigma) pair can be used only if it is one, otherwise it means that more than 1 Cluster Size
447  // distribution gave these values after being fitted with a Landau
448  TH2D* hFitR0CntAll;
449  fSimInput->GetObject("hClustSizeFitR0CntAll", hFitR0CntAll);
450  // Number of entries for each bin of the Sigma of a Landau Input distribution as function of
451  // the MPV and sigma values of a Landau fit on the corresponding Cluster size distribution
452  // The (R0, Input Sigma) pair can be used only if it is one, otherwise it means that more than 1 Cluster Size
453  // distribution gave these values after being fitted with a Landau
454  TH2D* hFitSigInCntAll;
455  fSimInput->GetObject("hClustSizeFitSigInCntAll", hFitSigInCntAll);
456  if (0 == hFitR0All || 0 == hFitSigInAll || 0 == hFitR0CntAll
457  || 0 == hFitSigInCntAll) {
458  LOG(error) << "CbmTofDigiBdfPar::GetLandauParFromBeamDataFit => Failed to "
459  "load Simulated Landau distribution values "
460  << " => Use default values from ASCII parameter file! Pointers: "
461  << hFitR0All << " " << hFitSigInAll << " " << hFitR0CntAll << " "
462  << hFitSigInCntAll << " ";
463  return kFALSE;
464  } // if( 0 == hFitR0All || 0 == hFitSigInAll || 0 == hFitR0CntAll || 0 == hFitSigInCntAll )
465 
466  // Prepare the TArrayD
469  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++) {
470  TFitResultPtr pResult = fh1ClusterSize[iSmType]->Fit("landau", "SN");
471  if (kTRUE == pResult->IsValid()) {
472  // Fit succedeed
473  Double_t dMpv = pResult->Parameter(1);
474  Double_t dSigma = pResult->Parameter(2);
475 
476  if (1 == hFitR0CntAll->GetBinContent(hFitR0CntAll->FindBin(dMpv, dSigma))
477  && 1
478  == hFitSigInCntAll->GetBinContent(
479  hFitSigInCntAll->FindBin(dMpv, dSigma))) {
480  // Unique Cluster size distribution giving these values
481  fdLandauMpv[iSmType] =
482  hFitR0All->GetBinContent(hFitR0All->FindBin(dMpv, dSigma));
483  fdLandauSigma[iSmType] =
484  hFitSigInAll->GetBinContent(hFitSigInAll->FindBin(dMpv, dSigma));
485  } // if unique pair
486  else {
487  // Pair is not unique
488  // => Do nothing, as default parameters should have been previously loaded
489  LOG(WARNING) << "CbmTofDigiBdfPar::GetLandauParFromBeamDataFit => Fit "
490  "values matching a non-unique param pair for Sm Type "
491  << iSmType << ": "
492  << hFitR0CntAll->GetBinContent(
493  hFitR0CntAll->FindBin(dMpv, dSigma));
494  LOG(WARNING) << " => Use default values from ASCII parameter file";
495  } // Nor unique pair
496  } // if( kTRUE == pResult->IsValid() )
497  else {
498  // Fit failed
499  // => Do nothing, as default parameters should have been previously loaded
500  LOG(WARNING) << "CbmTofDigiBdfPar::GetLandauParFromBeamDataFit => Fit "
501  "failed for Sm Type "
502  << iSmType
503  << " => Use default values from ASCII parameter file";
504  } // else of if( kTRUE == pResult->IsValid() )
505  } // for( Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType ++)
506 
507  gDirectory->cd(oldir->GetPath());
508 
509  return kTRUE;
510 }
511 /************************************************************************************/
512 // Geometry variables
513 Int_t CbmTofDigiBdfPar::GetNbSm(Int_t iSmType) const {
514  if (iSmType < fiNbSmTypes)
515  return fiNbSm[iSmType];
516  else
517  return 0;
518 }
519 Int_t CbmTofDigiBdfPar::GetNbRpc(Int_t iSmType) const {
520  if (iSmType < fiNbSmTypes)
521  return fiNbRpc[iSmType];
522  else
523  return 0;
524 }
525 Int_t CbmTofDigiBdfPar::GetNbGaps(Int_t iSmType, Int_t iRpc) const {
526  if (iSmType < fiNbSmTypes) {
527  if (iRpc < fiNbRpc[iSmType])
528  return fiNbGaps[iSmType][iRpc];
529  else
530  return 0;
531  } // if( iSmType < fiNbSmTypes )
532  else
533  return 0;
534 }
535 Double_t CbmTofDigiBdfPar::GetGapSize(Int_t iSmType, Int_t iRpc) const {
536  if (iSmType < fiNbSmTypes) {
537  if (iRpc < fiNbRpc[iSmType])
538  return fdGapSize[iSmType][iRpc];
539  else
540  return 0.0;
541  } // if( iSmType < fiNbSmTypes )
542  else
543  return 0.0;
544 }
545 Double_t
546 CbmTofDigiBdfPar::GetSigVel(Int_t iSmType, Int_t iSm, Int_t iRpc) const {
547  if (iSmType < fiNbSmTypes) {
548  if (iSm < fiNbSm[iSmType] && iRpc < fiNbRpc[iSmType])
549  return fdSigVel[iSmType][iSm * fiNbRpc[iSmType] + iRpc];
550  else
551  return 0.0;
552  } // if( iSmType < fiNbSmTypes )
553  else
554  return 0.0;
555 }
556 void CbmTofDigiBdfPar::SetSigVel(Int_t iSmType,
557  Int_t iSm,
558  Int_t iRpc,
559  Double_t vel) {
560  if (iSmType > -1 && iSmType < fiNbSmTypes) {
561  if (iSm < fiNbSm[iSmType] && iRpc < fiNbRpc[iSmType])
562  fdSigVel[iSmType][iSm * fiNbRpc[iSmType] + iRpc] = vel;
563  else {
564  LOG(error) << "Invalid iSm " << iSm << ", iRpc " << iRpc;
565  }
566  } else {
567  LOG(error) << "Invalid SMType " << iSmType;
568  }
569 }
570 Int_t CbmTofDigiBdfPar::GetNbChan(Int_t iSmType, Int_t iRpc) const {
571  if (iSmType < fiNbSmTypes) {
572  if (iRpc < fiNbRpc[iSmType])
573  return fiNbCh[iSmType][iRpc];
574  else
575  return 0;
576  } // if( iSmType < fiNbSmTypes )
577  else
578  return 0;
579 }
580 Int_t CbmTofDigiBdfPar::GetChanType(Int_t iSmType, Int_t iRpc) const {
581  if (iSmType < fiNbSmTypes) {
582  if (iRpc < fiNbRpc[iSmType])
583  return fiChType[iSmType][iRpc];
584  else
585  return 0;
586  } // if( iSmType < fiNbSmTypes )
587  else
588  return 0;
589 }
590 Int_t CbmTofDigiBdfPar::GetChanOrient(Int_t iSmType, Int_t iRpc) const {
591  if (iSmType < fiNbSmTypes) {
592  if (iRpc < fiNbRpc[iSmType])
593  return fiChOrientation[iSmType][iRpc];
594  else
595  return 0;
596  } // if( iSmType < fiNbSmTypes )
597  else
598  return 0;
599 }
600 // Beamtime variables
601 Int_t CbmTofDigiBdfPar::GetTypeInputMap(Int_t iSmType) const {
602  if (iSmType < fiNbSmTypes)
603  return fiSmTypeInpMapp[iSmType];
604  else
605  return 0;
606 }
607 Double_t CbmTofDigiBdfPar::GetEfficiency(Int_t iSmType) const {
608  if (iSmType < fiNbSmTypes)
609  return fdEfficiency[iSmType];
610  else
611  return 0;
612 }
613 Double_t CbmTofDigiBdfPar::GetGapEfficiency(Int_t iSmType, Int_t iRpc) const {
614  if (iSmType < fiNbSmTypes) {
615  if (iRpc < fiNbRpc[iSmType])
616  return fdGapsEfficiency[iSmType][iRpc];
617  else
618  return 0;
619  } // if( iSmType < fiNbSmTypes )
620  else
621  return 0;
622 }
623 Double_t CbmTofDigiBdfPar::GetResolution(Int_t iSmType) const {
624  if (iSmType < fiNbSmTypes)
625  return fdTimeResolution[iSmType];
626  else
627  return 0;
628 }
629 Double_t CbmTofDigiBdfPar::GetSystemResolution(Int_t iSmType) const {
630  if (iSmType < fiNbSmTypes)
631  return TMath::Sqrt(fdTimeResElec * fdTimeResElec
633  + fdTimeResolution[iSmType] * fdTimeResolution[iSmType]);
634  else
635  return 0;
636 }
637 TH1* CbmTofDigiBdfPar::GetClustSizeHist(Int_t iSmType) const {
638  if (iSmType < fiNbSmTypes)
639  return fh1ClusterSize[iSmType];
640  else
641  return 0;
642 }
643 TH1* CbmTofDigiBdfPar::GetClustTotHist(Int_t iSmType) const {
644  if (iSmType < fiNbSmTypes)
645  return fh1ClusterTot[iSmType];
646  else
647  return 0;
648 }
649 Double_t CbmTofDigiBdfPar::GetLandauMpv(Int_t iSmType) const {
650  if (iSmType < fiNbSmTypes)
651  return fdLandauMpv[iSmType];
652  else
653  return 0;
654 }
655 Double_t CbmTofDigiBdfPar::GetLandauSigma(Int_t iSmType) const {
656  if (iSmType < fiNbSmTypes)
657  return fdLandauSigma[iSmType];
658  else
659  return 0;
660 }
662  if (-1 < fdMaxTimeDistClust)
663  return fdMaxTimeDistClust;
664  else
665  return 5 * 0.08;
666 }
668  // LOG(info)<<" _ _ _ _ _ _ _ _ _ _ ";
669  // LOG(info)<<" / \_/ \_/ \_/ \_/ \_/ \_/ \_/ \_/ \_/ \ ";
670  // LOG(info)<<" \_/ \_/ \_/ \_/ \_/ \_/ \_/ \_/ \_/ \_/ ";
671  LOG(info) << " _ _ _ _ _ _ _ _ _ _ ";
672  LOG(info) << " |_|+|_|+|_|+|_|+|_|+|_|+|_|+|_|+|_|+|_| ";
673  LOG(info) << "Parameter values in CbmTofDigiBdfPar: ";
674  LOG(info) << "=> MC data usage: ";
675  if (kTRUE == fbUseOnlyPrim)
676  LOG(info) << " Tracks used: Only Primary!";
677  else
678  LOG(info)
679  << " Tracks used: All (Primary + Secondary)";
680  // Fee properties
681  LOG(info) << "=> FEE properties: ";
682  LOG(info) << " FEE gain standard deviation: "
683  << 100.0 * fdFeeGainSigma << " [%]";
684  LOG(info) << " FEE Threshold on ToT: " << fdFeeTotThr
685  << "[ns]";
686  LOG(info) << " FEE channel time resolution: " << fdTimeResElec
687  << "[ns]";
688  LOG(info) << " Start channel time resolution: " << fdTimeResStart
689  << "[ns]";
690  LOG(info) << " FEE deadtime: " << fdDeadtime << "[ns]";
691 
692  // Geometry variables
693  LOG(info) << "=> Geometry variables: ";
694  LOG(info) << " Signal propa. speed on readout el.: " << fdSignalPropSpeed
695  << " [cm/ns]";
696  LOG(info) << " Number of Super Modules (SM) Types: " << fiNbSmTypes;
697 
698  TString sIndex = " Sm Type |-- ";
699  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++)
700  sIndex += Form("%3d ", iSmType);
701  LOG(info) << sIndex;
702 
703  TString sSmNb = " Nb SM for each SM type: |-> ";
704  TString sRpcNb = " Nb Rpc for each SM type: |-> ";
705  sIndex = " Rpc index |-- ";
706  Int_t iMaxRpcNb = 0;
707 
708  TString* sGapsNb = new TString[fiNbSmTypes];
709  TString* sGapsSz = new TString[fiNbSmTypes];
710  TString* sSigVel = new TString[fiNbSmTypes];
711  TString* sChNb = new TString[fiNbSmTypes];
712  TString* sChType = new TString[fiNbSmTypes];
713  TString* sChOrient = new TString[fiNbSmTypes];
714 
715  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++) {
716  sSmNb += Form("%3d ", GetNbSm(iSmType));
717  sRpcNb += Form("%3d ", GetNbRpc(iSmType));
718  sGapsNb[iSmType] = Form(" Nb of Gaps in SM type %3d:|-> ", iSmType);
719  sGapsSz[iSmType] = Form(" Gap Size(mm) in SM type %3d:|-> ", iSmType);
720  sSigVel[iSmType] = Form(" SigVel(cm/ps) in SM type %3d:|-> ", iSmType);
721  sChNb[iSmType] = Form(" Nb of Chan in SM type %3d:|-> ", iSmType);
722  sChType[iSmType] = Form(" Chan Type in SM type %3d:|-> ", iSmType);
723  sChOrient[iSmType] = Form(" Chan orient. in SM type %3d:|-> ", iSmType);
724  if (iMaxRpcNb < fiNbRpc[iSmType]) iMaxRpcNb = fiNbRpc[iSmType];
725 
726 
727  for (Int_t iRpc = 0; iRpc < fiNbRpc[iSmType]; iRpc++) {
728  sGapsNb[iSmType] += Form("%3d ", GetNbGaps(iSmType, iRpc));
729  sGapsSz[iSmType] += Form("%3.2f ", GetGapSize(iSmType, iRpc));
730  for (Int_t iSm = 0; iSm < fiNbSm[iSmType]; iSm++)
731  sSigVel[iSmType] += Form("%4.3f ", GetSigVel(iSmType, iSm, iRpc));
732  sChNb[iSmType] += Form("%3d ", GetNbChan(iSmType, iRpc));
733  if (1 == GetChanType(iSmType, iRpc))
734  sChType[iSmType] += "pad ";
735  else
736  sChType[iSmType] += "str ";
737  if (1 == GetChanOrient(iSmType, iRpc))
738  sChOrient[iSmType] += "hor ";
739  else
740  sChOrient[iSmType] += "ver ";
741  }
742  } // for( Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType ++)
743  for (Int_t iRpc = 0; iRpc < iMaxRpcNb; iRpc++)
744  sIndex += Form("%3d ", iRpc);
745 
746  LOG(info) << sSmNb;
747  LOG(info) << sRpcNb;
748  LOG(info) << sIndex;
749  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++) {
750  LOG(info) << sGapsNb[iSmType];
751  LOG(info) << sGapsSz[iSmType];
752  LOG(debug) << sSigVel[iSmType];
753  LOG(info) << sChNb[iSmType];
754  LOG(info) << sChType[iSmType];
755  LOG(info) << sChOrient[iSmType];
756  } // for( Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType ++)
757 
758  // Beamtime variables
759  LOG(info) << "=> Beamtime variables: ";
760  LOG(info) << " Beamtime data input file name: " << fsBeamInputFile;
761  LOG(info) << " Beamtime data calib file name: " << fsBeamCalibFile;
762  switch (fiClusterRadiusModel) {
763  case -1:
764  LOG(info)
765  << " Cluster radius model: Single value at 0.0002 cm";
766  break;
767  case 0:
768  LOG(info) << " Cluster radius model: Single value using "
769  "beamtime mean cluster size";
770  break;
771  case 1:
772  LOG(info) << " Cluster radius model: Landau Distrib. "
773  "using values from parameter file";
774  break;
775  case 2:
776  LOG(info) << " Cluster radius model: Landau Distrib. "
777  "using values from fit on beam data";
778  break;
779  default:
780  LOG(info)
781  << " Cluster radius model: None, this should not be!!";
782  break;
783  } // switch( fiClusterRadiusModel )
784  sIndex = " Sm Type |-- ";
785  TString sTypeMap = " SM type in file: |-> ";
786  TString sEffType = " Efficiency: |-> ";
787  TString sResType = " Resolution: |-> ";
788  TString sSysResTy = " => Sys. Res.: |-> ";
789  TString sLandMpv = " Landau MPV: |-> ";
790  TString sLandSig = " Landau Sigma: |-> ";
791  for (Int_t iSmType = 0; iSmType < fiNbSmTypes; iSmType++) {
792  sIndex += Form("%3d ", iSmType);
793  sTypeMap += Form("%3d ", GetTypeInputMap(iSmType));
794  sEffType += Form("%1.2f ", fdEfficiency[iSmType]);
795  sResType += Form("%2.3f ", fdTimeResolution[iSmType]);
796  sSysResTy += Form("%2.3f ", GetSystemResolution(iSmType));
797  switch (fiClusterRadiusModel) {
798  case 1:
799  case 2:
800  sLandMpv += Form("%2.3f ", fdLandauMpv[iSmType]);
801  sLandSig += Form("%2.3f ", fdLandauSigma[iSmType]);
802  break;
803  default: break;
804  } // switch( fiClusterRadiusModel )
805  }
806  LOG(info) << sIndex;
807  LOG(info) << sTypeMap;
808  LOG(info) << sEffType;
809  LOG(info) << sResType;
810  LOG(info) << sSysResTy;
811  switch (fiClusterRadiusModel) {
812  case 1:
813  case 2:
814  LOG(info) << sLandMpv;
815  LOG(info) << sLandSig;
816  break;
817  default: break;
818  } // switch( fiClusterRadiusModel )
819 
820  LOG(info) << "=> Digitization variables: ";
821  // Digi type to use
822  if (kTRUE == fbUseExpDigi)
823  LOG(info)
824  << " Output ToF digi type: Expanded (Double values)";
825  else
826  LOG(info) << " Output ToF digi type: Re-digitized (Encoding "
827  "in 64 bits)";
828  if (kTRUE == fbOneGapTrack)
829  LOG(info) << " Avoid multi gaps digi for same trk: ON";
830  else
831  LOG(info) << " Avoid multi gaps digi for same trk: OFF";
832 
833  // Cluster Model to use
834  switch (fiClusterModel) {
835  case 0:
836  LOG(info) << " Cluster model: Use directly cluster "
837  "size (no cluster radius)";
838  break;
839  case 1:
840  LOG(info) << " Cluster model: Flat disc as charge "
841  "distribution + geometrical overlap";
842  break;
843  case 2:
844  LOG(info) << " Cluster model: 2D gaussian as "
845  "charge distribution + integral";
846  break;
847  default:
848  LOG(info)
849  << " Cluster model: None, this should not be!!";
850  break;
851  } // switch( fiClusterModel )
852  LOG(info) << "=> Simple clusterizer parameters: ";
853  if (kTRUE == fbMulUseTrackId)
854  LOG(info) << " Variable used for multiplicity cnt: Track ID";
855  else
856  LOG(info) << " Variable used for multiplicity cnt: TofPoint ID";
857  if (-1 < fdMaxTimeDistClust)
858  LOG(info) << " Maximal time dist. to last chan.: " << GetMaxTimeDist()
859  << " [ns]";
860  else
861  LOG(info) << " Maximal time dist. to last chan.: Use 5*Nom. Syst. Res. "
862  "(80 ps) = "
863  << GetMaxTimeDist() << " [ns]";
864 
865  LOG(info) << " Maximal dist. along ch to last one: " << fdMaxSpaceDistClust
866  << " [cm]";
867 
868 
869  delete[] sGapsNb;
870  delete[] sGapsSz;
871  delete[] sChNb;
872  delete[] sChType;
873  delete[] sChOrient;
874 }
875 
876 Int_t CbmTofDigiBdfPar::GetNbDet() const { return fiDetUId.GetSize(); }
877 
878 Int_t CbmTofDigiBdfPar::GetDetInd(Int_t iAddr) {
879 
880  const Int_t DetMask = 4194303;
881  /*
882  Int_t iInd=0;
883  while (fiDetUId[iInd] != (iAddr & DetMask) && iInd<fiDetUId.GetSize()) iInd++; //FIXME, inefficient and using a numerical constant!
884  if(iInd == fiDetUId.GetSize()){
885  LOG(error)<<Form(" detector id 0x%08x unknown ",iAddr) ;
886  iInd=0;
887  }
888  return iInd;
889  */
890  return fMapDetInd[(iAddr & DetMask)];
891 }
892 
893 Int_t CbmTofDigiBdfPar::GetDetUId(Int_t iInd) { return fiDetUId[iInd]; }
CbmTofDigiBdfPar.h
CbmTofDigiBdfPar::fdLandauMpv
TArrayD fdLandauMpv
Definition: CbmTofDigiBdfPar.h:154
CbmTofDigiBdfPar::fiNbGaps
std::vector< TArrayI > fiNbGaps
Definition: CbmTofDigiBdfPar.h:132
CbmTofDigiBdfPar::fdDeadtime
Double_t fdDeadtime
Definition: CbmTofDigiBdfPar.h:125
CbmTofDigiBdfPar::printParams
void printParams()
Definition: CbmTofDigiBdfPar.cxx:667
CbmTofDigiBdfPar::fbMulUseTrackId
Bool_t fbMulUseTrackId
Definition: CbmTofDigiBdfPar.h:158
CbmTofDigiBdfPar::fbUseOnlyPrim
Bool_t fbUseOnlyPrim
Definition: CbmTofDigiBdfPar.h:112
CbmTofDigiBdfPar::fdGapSize
std::vector< TArrayD > fdGapSize
Definition: CbmTofDigiBdfPar.h:133
CbmTofDigiBdfPar::fiDetUId
TArrayI fiDetUId
Definition: CbmTofDigiBdfPar.h:140
CbmTofDigiBdfPar::GetClustSizeHist
TH1 * GetClustSizeHist(Int_t iSmType) const
Definition: CbmTofDigiBdfPar.cxx:637
CbmTofDigiBdfPar::fiNbRpc
TArrayI fiNbRpc
Definition: CbmTofDigiBdfPar.h:131
CbmTofDigiBdfPar::fh1ClusterTot
std::vector< TH1 * > fh1ClusterTot
Definition: CbmTofDigiBdfPar.h:153
CbmTofDigiBdfPar::GetResolution
Double_t GetResolution(Int_t iSmType) const
Definition: CbmTofDigiBdfPar.cxx:623
CbmTofDigiBdfPar::fMapDetInd
std::map< Int_t, Int_t > fMapDetInd
Definition: CbmTofDigiBdfPar.h:141
CbmTofDigiBdfPar::GetNbChan
Int_t GetNbChan(Int_t iSmType, Int_t iRpc) const
Definition: CbmTofDigiBdfPar.cxx:570
ClassImp
ClassImp(CbmTofDigiBdfPar) CbmTofDigiBdfPar
Definition: CbmTofDigiBdfPar.cxx:24
CbmTofDigiBdfPar::GetGapEfficiency
Double_t GetGapEfficiency(Int_t iSmType, Int_t iRpc) const
Definition: CbmTofDigiBdfPar.cxx:613
CbmTofDigiBdfPar::fiNbCh
std::vector< TArrayI > fiNbCh
Definition: CbmTofDigiBdfPar.h:136
i
int i
Definition: L1/vectors/P4_F32vec4.h:25
CbmTofDigiBdfPar::GetNbRpc
Int_t GetNbRpc(Int_t iSmType) const
Definition: CbmTofDigiBdfPar.cxx:519
CbmTofDigiBdfPar::fiNbSm
TArrayI fiNbSm
Definition: CbmTofDigiBdfPar.h:130
CbmTofDigiBdfPar::GetChanOrient
Int_t GetChanOrient(Int_t iSmType, Int_t iRpc) const
Definition: CbmTofDigiBdfPar.cxx:590
CbmTofDigiBdfPar::GetSigVel
Double_t GetSigVel(Int_t iSmType, Int_t iSm, Int_t iRpc) const
Definition: CbmTofDigiBdfPar.cxx:546
CbmTofDigiBdfPar::GetDetInd
Int_t GetDetInd(Int_t iAddr)
Definition: CbmTofDigiBdfPar.cxx:878
CbmTofDigiBdfPar::fdFeeGainSigma
Double_t fdFeeGainSigma
Definition: CbmTofDigiBdfPar.h:121
CbmTofDigiBdfPar::fdSignalPropSpeed
Double_t fdSignalPropSpeed
Definition: CbmTofDigiBdfPar.h:128
CbmTofDigiBdfPar::putParams
void putParams(FairParamList *)
Definition: CbmTofDigiBdfPar.cxx:88
CbmTofDigiBdfPar::fbUseExpDigi
Bool_t fbUseExpDigi
Definition: CbmTofDigiBdfPar.h:110
CbmTofDigiBdfPar::fiClusterRadiusModel
Int_t fiClusterRadiusModel
Definition: CbmTofDigiBdfPar.h:147
CbmTofDigiBdfPar::fdTimeResElec
Double_t fdTimeResElec
Definition: CbmTofDigiBdfPar.h:123
CbmTofDigiBdfPar::getParams
Bool_t getParams(FairParamList *)
Definition: CbmTofDigiBdfPar.cxx:131
CbmTofDigiBdfPar::fiChType
std::vector< TArrayI > fiChType
Definition: CbmTofDigiBdfPar.h:137
CbmTofDigiBdfPar::GetTypeInputMap
Int_t GetTypeInputMap(Int_t iSmType) const
Definition: CbmTofDigiBdfPar.cxx:601
CbmTofDigiBdfPar::fdFeeTotThr
Double_t fdFeeTotThr
Definition: CbmTofDigiBdfPar.h:122
CbmTofAddress.h
CbmTofDigiBdfPar::GetChanType
Int_t GetChanType(Int_t iSmType, Int_t iRpc) const
Definition: CbmTofDigiBdfPar.cxx:580
CbmTofDigiBdfPar
Parameters class for the CBM ToF digitizer using beam data distributions.
Definition: CbmTofDigiBdfPar.h:30
CbmTofDigiBdfPar::GetLandauSigma
Double_t GetLandauSigma(Int_t iSmType) const
Definition: CbmTofDigiBdfPar.cxx:655
CbmTofDigiBdfPar::fdSigVel
std::vector< TArrayD > fdSigVel
Definition: CbmTofDigiBdfPar.h:135
CbmTofDigiBdfPar::~CbmTofDigiBdfPar
~CbmTofDigiBdfPar(void)
Definition: CbmTofDigiBdfPar.cxx:68
CbmTofDigiBdfPar::fbOneGapTrack
Bool_t fbOneGapTrack
Definition: CbmTofDigiBdfPar.h:115
CbmTofDigiBdfPar::GetLandauMpv
Double_t GetLandauMpv(Int_t iSmType) const
Definition: CbmTofDigiBdfPar.cxx:649
CbmTofDigiBdfPar::GetLandauParFromBeamDataFit
Bool_t GetLandauParFromBeamDataFit()
Definition: CbmTofDigiBdfPar.cxx:424
CbmTofDigiBdfPar::GetGapSize
Double_t GetGapSize(Int_t iSmType, Int_t iRpc) const
Definition: CbmTofDigiBdfPar.cxx:535
CbmTofDigiBdfPar::fiSmTypeInpMapp
TArrayI fiSmTypeInpMapp
Definition: CbmTofDigiBdfPar.h:148
CbmTofDigiBdfPar::fdGapsEfficiency
std::vector< TArrayD > fdGapsEfficiency
Definition: CbmTofDigiBdfPar.h:150
CbmTofDigiBdfPar::GetSystemResolution
Double_t GetSystemResolution(Int_t iSmType) const
Definition: CbmTofDigiBdfPar.cxx:629
CbmTofDigiBdfPar::fdTimeResolution
TArrayD fdTimeResolution
Definition: CbmTofDigiBdfPar.h:151
CbmTofDigiBdfPar::fdEfficiency
TArrayD fdEfficiency
Definition: CbmTofDigiBdfPar.h:149
CbmTofDigiBdfPar::GetNbDet
Int_t GetNbDet() const
Definition: CbmTofDigiBdfPar.cxx:876
CbmTofDigiBdfPar::fiClusterModel
Int_t fiClusterModel
Definition: CbmTofDigiBdfPar.h:118
CbmTofDigiBdfPar::fiNbSmTypes
Int_t fiNbSmTypes
Definition: CbmTofDigiBdfPar.h:129
CbmTofDigiBdfPar::fdMaxTimeDistClust
Double_t fdMaxTimeDistClust
Definition: CbmTofDigiBdfPar.h:159
CbmTofDigiBdfPar::fsBeamInputFile
TString fsBeamInputFile
Definition: CbmTofDigiBdfPar.h:144
CbmTofDigiBdfPar::fsBeamCalibFile
TString fsBeamCalibFile
Definition: CbmTofDigiBdfPar.h:145
CbmTofDigiBdfPar::SetSigVel
void SetSigVel(Int_t iSmType, Int_t iSm, Int_t iRpc, Double_t dvel)
Definition: CbmTofDigiBdfPar.cxx:556
CbmTofDigiBdfPar::fdTimeResStart
Double_t fdTimeResStart
Definition: CbmTofDigiBdfPar.h:124
CbmTofAddress::GetUniqueAddress
static UInt_t GetUniqueAddress(UInt_t Sm, UInt_t Rpc, UInt_t Channel, UInt_t Side=0, UInt_t SmType=0)
Definition: CbmTofAddress.h:124
CbmTofDigiBdfPar::GetDetUId
Int_t GetDetUId(Int_t iDet)
Definition: CbmTofDigiBdfPar.cxx:893
CbmTofDigiBdfPar::fiChOrientation
std::vector< TArrayI > fiChOrientation
Definition: CbmTofDigiBdfPar.h:139
CbmTofDigiBdfPar::clear
void clear(void)
Definition: CbmTofDigiBdfPar.cxx:74
CbmTofDigiBdfPar::fdLandauSigma
TArrayD fdLandauSigma
Definition: CbmTofDigiBdfPar.h:155
CbmTofDigiBdfPar::GetMaxTimeDist
Double_t GetMaxTimeDist() const
Definition: CbmTofDigiBdfPar.cxx:661
CbmTofDigiBdfPar::fdMaxSpaceDistClust
Double_t fdMaxSpaceDistClust
Definition: CbmTofDigiBdfPar.h:160
CbmTofDigiBdfPar::ClearHistos
void ClearHistos()
Definition: CbmTofDigiBdfPar.cxx:79
CbmTofDigiBdfPar::GetEfficiency
Double_t GetEfficiency(Int_t iSmType) const
Definition: CbmTofDigiBdfPar.cxx:607
CbmTofDigiBdfPar::GetClustTotHist
TH1 * GetClustTotHist(Int_t iSmType) const
Definition: CbmTofDigiBdfPar.cxx:643
CbmTofDigiBdfPar::CbmTofDigiBdfPar
CbmTofDigiBdfPar(const char *name="CbmTofDigiBdfPar", const char *title="BDF Digitization parameters for the TOF detector", const char *context="TestDefaultContext")
DetMask
const Int_t DetMask
Definition: CbmTofAnaTestbeam.cxx:75
CbmTofDigiBdfPar::fh1ClusterSize
std::vector< TH1 * > fh1ClusterSize
Definition: CbmTofDigiBdfPar.h:152
CbmTofDigiBdfPar::GetNbGaps
Int_t GetNbGaps(Int_t iSmType, Int_t iRpc) const
Definition: CbmTofDigiBdfPar.cxx:525
CbmTofDigiBdfPar::GetNbSm
Int_t GetNbSm(Int_t iSmType) const
Definition: CbmTofDigiBdfPar.cxx:513
CbmTofDigiBdfPar::LoadBeamtimeHistos
Bool_t LoadBeamtimeHistos()
Definition: CbmTofDigiBdfPar.cxx:314