CbmRoot
CbmFastSim.cxx
Go to the documentation of this file.
1 //
3 // FastSim
4 //
7 
8 //C++ class headers
9 #include <fstream>
10 #include <string>
11 
12 //Fair class headers
13 #include "FairLogger.h"
14 #include "FairRootManager.h"
15 //#include "FairRunAna.h"
16 #include "FairRun.h"
17 #include "FairRuntimeDb.h"
18 
19 //CBM class headers
20 #include "CbmMCTrack.h"
21 #include "CbmStack.h"
22 #include "PairAnalysisHelper.h"
23 
24 //ROOT class headers
25 #include <TClonesArray.h>
26 #include <TDatabasePDG.h>
27 #include <TF1.h>
28 #include <TH1.h>
29 #include <TH2.h>
30 #include <TH3.h>
31 #include <THnBase.h>
32 #include <THnSparse.h>
33 #include <TMath.h>
34 #include <TMatrixD.h>
35 #include <TPDGCode.h>
36 #include <TParticle.h>
37 #include <TParticlePDG.h>
38 #include <TProfile3D.h>
39 #include <TRandom3.h>
40 #include <TVector3.h>
41 #include <TVectorD.h>
42 #include <TVirtualMC.h>
43 
44 
45 #include "CbmFastSim.h"
46 
47 
48 using std::cout;
49 using std::endl;
50 using std::ifstream;
51 using std::string;
52 
53 
54 // ----- Default constructor -------------------------------------------
56  : FairTask("Cbm Fast Simulation")
57  , fRand(new TRandom3(0))
58  , fdbPdg(TDatabasePDG::Instance()) {
62  for (Int_t idx = 0; idx < fgkNParts; idx++) {
63  fEFF[idx] = NULL;
64  fEFFdenom[idx] = NULL;
65  for (Int_t j = 0; j < 10; j++)
66  fEFFproj[idx][j] = NULL;
67 
68  fRFp[idx] = NULL;
69  }
70 }
71 // -------------------------------------------------------------------------
72 
73 // ----- Destructor ----------------------------------------------------
78  FairRootManager* fManager = FairRootManager::Instance();
79  fManager->Write();
80 
81  if (fFastTracks) {
82  fFastTracks->Delete();
83  delete fFastTracks;
84  }
85  if (fRand) delete fRand;
86 }
87 // -------------------------------------------------------------------------
88 
89 
94  fFastTracks = new TClonesArray("TParticle");
95  FairRootManager::Instance()->Register(
96  "FastTrack", "FastSim", fFastTracks, kTRUE);
97 }
98 
99 
100 // ----- Public method Init --------------------------------------------
102 
103 // ----- Public method Init --------------------------------------------
104 InitStatus CbmFastSim::Init() {
105 
107  for (Int_t idx = 0; idx < fgkNParts; idx++) {
108 
109  if (fEFF[idx] && fEFFdenom[idx]) {
110  Info("Init",
111  "Prepare efficiency map and projections for index %d and method %d",
112  idx,
113  fMethod);
114  // fEFF[idx]->Print("a");
115 
116  // first get projections and calculate integrated effiencies
117  // NOTE: if statistic are bad these are used
118  for (Int_t j = 0; j < fEFF[idx]->GetNdimensions(); j++) {
119 
121  switch (fMethod) {
122  case kIgnoreFluct:
124  continue;
125  break;
126  case kLastDim:
128  if (j < (fEFF[idx]->GetNdimensions() - 1)) continue;
129  case kAverage:
131  case kFactorize: {
133  fEFFproj[idx][j] = fEFF[idx]->Projection(j, "E");
134  fEFFproj[idx][j]->SetName(Form("NOMidx%dj%d", idx, j));
135 
136  TH1* den = fEFFdenom[idx]->Projection(j, "E");
137  den->SetName(Form("DENidx%dj%d", idx, j));
138  fEFFproj[idx][j]->Divide(den);
139 
140  // factoring using last dimension!=1 histograms
141  if (fMethod == kFactorize && j > 0)
142  fEFFproj[idx][j]->Scale(1.
143  / fEFFproj[idx][j]->GetBinContent(
144  fEFFproj[idx][j]->GetMaximumBin()));
145 
146  if (den) delete den;
147  } break;
148  case kInterpolate: {
150  if (j > 0) continue;
151  fEFFproj[idx][j] = fEFF[idx]->Projection(j, 1, 2, "E");
152  fEFFproj[idx][j]->SetName(Form("NOMidx%dj%d12", idx, j));
153 
154  TH1* den = fEFFdenom[idx]->Projection(j, 1, 2, "E");
155  den->SetName(Form("DENidx%dj%d12", idx, j));
156  fEFFproj[idx][j]->Divide(den);
157 
158  if (den) delete den;
159  } break;
160  default: break;
161  }
162  }
163 
164  // calculate main effiency map
165  fEFF[idx]->Divide(fEFFdenom[idx]);
166 
167  } // endif: nominator
168 
169 
171  if (fRFp[idx]) {
172  Info("Init", "Momentum smearing map found for index %d", idx);
173  PairAnalysisHelper::CumulateSlicesX(fRFp[idx], kFALSE, kTRUE);
174  }
175 
176  } //end: particle species
177 
178  // Create and register output array
179  Register();
180 
181  return kSUCCESS;
182 }
183 
185 
186  // Get run and runtime database
187  FairRun* run = FairRun::Instance();
188  if (!run) Fatal("SetParContainers", "No analysis run");
189 
190  FairRuntimeDb* db = run->GetRuntimeDb();
191  if (!db) Fatal("SetParContainers", "No runtime database");
192 }
193 // -------------------------------------------------------------------------
194 
195 void CbmFastSim::SetSeed(unsigned int seed) {
199  fRand->SetSeed(seed);
200 }
201 
202 // ----- Public method Finish --------------------------------------------
207 }
208 // -------------------------------------------------------------------------
209 
210 // ----- Public method Exec --------------------------------------------
211 void CbmFastSim::Exec(Option_t* opt) {
215 
217  CbmStack* fStack = (CbmStack*) gMC->GetStack();
218  Int_t nTracks = fStack->GetNtrack();
219  // Info("Exec","number of tracks **** %d",Tracks);
220 
222  if (fFastTracks->GetEntriesFast() != 0) fFastTracks->Clear("C");
223 
225  TClonesArray& chrgCandidates = *fFastTracks; // Charged Candidates
226 
227  Int_t iKeep = 0;
228  Int_t iConv = 0;
229 
231  Double_t* vals = new Double_t[4]; //fEFF->GetNdimensions()];
232  Int_t* coord = new Int_t[4];
233 
235  for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
236 
238  TParticle* t = fStack->GetParticle(iTrack);
239  int pdg = t->GetPdgCode();
240  // t->Print();
241 
243  if (!fdbPdg->GetParticle(pdg)) continue;
244  if (TMath::Abs(fdbPdg->GetParticle(pdg)->Charge()) < 3) continue;
245 
248  vals[0] = t->Pt();
249  vals[1] = t->Theta();
250  vals[2] = TVector2::Phi_0_2pi(t->Phi());
251  vals[3] = t->Vz();
252 
255  if (TMath::Abs(pdg) == 11 && fStack->GetParticle(t->GetMother(0))
256  && fStack->GetParticle(t->GetMother(0))->GetPdgCode() == 22) {
257  if (!PassEfficiencyFilter(22, vals, coord)) continue;
258  iConv++;
259  } else {
260  if (!PassEfficiencyFilter(pdg, vals, coord)) continue;
261  }
262 
264  TLorentzVector p4(t->Px(), t->Py(), t->Pz(), t->Energy());
265 
267  // TVector3 stvtx(t->Vx(),t->Vy(),t->Vz());
268 
270  if (!Smearing(&p4, pdg)) { /*ERROR*/
271  continue;
272  }
273 
276  new ((*fFastTracks)[iKeep]) TParticle(pdg,
277  t->GetStatusCode(),
278  iTrack,
279  0,
280  0,
281  0,
282  p4.Px(),
283  p4.Py(),
284  p4.Pz(),
285  p4.Energy(),
286  t->Vx(),
287  t->Vy(),
288  t->Vz(),
289  0);
290 
291  iKeep++;
292 
293  } //trackloop
294 
295  // cleanup
296  delete[] vals;
297  delete[] coord;
298 
299  Info("Exec", "stack size: \t %d \n", fStack->GetNtrack());
300  Info("Exec", "kept tracks: \t %d \n", iKeep);
301  Info("Exec", "kept converison tracks: \t %d \n", iConv);
302 }
303 
304 Bool_t
305 CbmFastSim::PassEfficiencyFilter(Int_t pdg, Double_t* vals, Int_t* coord) {
306  //
307  // check if partilce passed the reconstruction cuts
308  //
309  Int_t idx = -1;
310  switch (TMath::Abs(pdg)) {
311  case kGamma: idx = CbmFastSim::kGam; break;
312  case kElectron: idx = CbmFastSim::kEle; break;
313  case kMuonMinus: idx = CbmFastSim::kMuo; break;
314  case kPiPlus: idx = CbmFastSim::kPio; break;
315  case kKPlus: idx = CbmFastSim::kKao; break;
316  case kProton: idx = CbmFastSim::kPro; break;
317  default: /*Warning("PassEfficiencyFilter","Pdg code %d not supported, particle will not be processed",pdg);*/
318  return kFALSE;
319  }
320 
322  if (!fEFF[idx]) {
323  Warning("PassEfficiencyFilter", "No filter for pdg code %d", pdg);
324  return kFALSE;
325  }
326 
328  // Info("PassEfficiencyFilter","Efficiency for pt:%.2e theta:%.2e phi:%.2e", vals[0],vals[1],vals[2]);
329 
331  Long_t bin = fEFF[idx]->GetBin(vals, kFALSE);
332 
334  Double_t eff = 0.0;
335  Double_t err = 0.0;
336  if (bin > 0) {
337  eff =
338  fEFF[idx]->GetBinContent(bin, coord); // write bin coordinates into coord
339 
341  if (fMethod != kIgnoreFluct) {
342 
343  err = fEFF[idx]->GetBinError(bin);
344  // if (eff>0. && err/eff>=0.60 ) return kFALSE;
345  if (eff > 0. && err / eff >= 0.10) {
346  Double_t effm = GetEfficiency(idx, vals, coord); // vals was NULL
347  // Info("PassEfficiencyFilter","statistical error %f too small (eff:%.2e, new:%.2e)",(eff>0.?err/eff:0.),eff,effm);
348  eff = effm;
349  }
350  }
351  }
352  // else {
353  // eff= GetEfficiency(idx, vals, NULL);
354  // Info("PassEfficiencyFilter","No entry in map, caculate via projections: %.2e (mean eff.)",eff);
355  // }
356 
357  Double_t rndm = fRand->Rndm();
358 
359  if (rndm > eff)
360  return kFALSE;
361  else
362  return kTRUE;
363 }
364 
365 Bool_t CbmFastSim::Smearing(TLorentzVector* p4, Int_t pdg) {
370 
371  Int_t idx = -1;
372  switch (TMath::Abs(pdg)) {
373  case kElectron: idx = CbmFastSim::kEle; break;
374  case kMuonMinus: idx = CbmFastSim::kMuo; break;
375  case kPiPlus: idx = CbmFastSim::kPio; break;
376  case kKPlus: idx = CbmFastSim::kKao; break;
377  case kProton: idx = CbmFastSim::kPro; break;
378  default: return kFALSE;
379  }
380 
382  Double_t tP = p4->P();
383  Double_t tTheta = p4->Theta();
384  Double_t tPhi = TVector2::Phi_0_2pi(p4->Phi());
385  Double_t mass = fdbPdg->GetParticle(pdg)->Mass();
386 
388  Double_t sP = AnalyzeMap(fRFp[idx], tP);
389  Double_t sTheta = tTheta;
390  Double_t sPhi = tPhi;
391 
393  Double_t sPx = sP * TMath::Sin(sTheta) * TMath::Cos(sPhi);
394  Double_t sPy = sP * TMath::Sin(sTheta) * TMath::Sin(sPhi);
395  Double_t sPz = sP * TMath::Cos(sTheta);
396 
397  Double_t sPt = TMath::Sqrt(sPx * sPx + sPy * sPy);
398  Double_t eta = -TMath::Log(TMath::Tan(sTheta / 2));
399 
400 
401  // p4->Print();
402  p4->SetPtEtaPhiM(sPt, eta, sPhi, mass);
403  // p4->Print();
404  // printf("\n\n");
405 
406  return kTRUE;
407 }
408 
409 Double_t CbmFastSim::AnalyzeMap(TH2F* map, Double_t refValue) {
414 
416  if (!map) return refValue;
417 
419  TAxis* xRec = map->GetXaxis();
420  TAxis* yGen = map->GetYaxis();
421 
423  Int_t binGen = yGen->FindBin(refValue);
424 
429  if (binGen < 1 || binGen > yGen->GetNbins()) {
430  return (fRand->Gaus(refValue, refValue * 0.02));
431  }
432 
434  Float_t* arrayh = map->GetArray();
435  Int_t offset = (yGen->GetNbins() + 2) * (binGen) + 1;
436 
438  Float_t r1 = fRand->Rndm();
439  Int_t binRec =
440  TMath::BinarySearch(xRec->GetNbins() + 2, arrayh + offset - 1, r1);
441 
442  Double_t smearedValue = xRec->GetBinLowEdge(binRec + 1);
443  if (r1 > arrayh[offset - 1 + binRec])
444  smearedValue +=
445  xRec->GetBinWidth(binRec + 1) * (r1 - arrayh[offset - 1 + binRec])
446  / (arrayh[offset - 1 + binRec + 1] - arrayh[offset - 1 + binRec]);
447 
448  // printf("binGen: %d, binRec: %d , offset: %d , \t r1: %f \t in: %f --> out: %f \n",
449  // binGen,binRec,offset, r1,refValue,smearedValue);
450 
451  return smearedValue;
452 }
453 
454 
456  THnBase* denom,
457  EParticleType part) {
461 
462  if (!nom || !denom)
463  Fatal("SetLookupEfficiency", "Either nominator or denominator is NULL ");
464  fEFF[part] = nom;
465  fEFFdenom[part] = denom;
466  /*
467  // first get projections and calculate integrated effiencies
468  // NOTE: if statistic are bad these are used
469  for (Int_t j = 0; j < nom->GetNdimensions(); ++j) {
470 
472  switch(fMethod) {
473  case kIgnoreFluct:
475  continue;
476  break;
477  case kLastDim:
479  if( j<(nom->GetNdimensions()-1) ) continue;
480  break;
481  case kAverage:
483  case kFactorize: {
485  fEFFproj[part][j] = nom->Projection(j,"E");
486  fEFFproj[part][j]->SetName(Form("NOMidx%dj%d",part,j));
487 
488  TH1 *den=denom->Projection(j,"E");
489  den->SetName(Form("DENidx%dj%d",part,j));
490  fEFFproj[part][j]->Divide( den );
491 
492  // factoring using last dimension!=1 histograms
493  if(fMethod==kFactorize && j>0) fEFFproj[part][j]->Scale( 1./fEFFproj[part][j]->GetBinContent( fEFFproj[part][j]->GetMaximumBin() ) );
494 
495  if(den) delete den;
496  }
497  break;
498  case kInterpolate: {
500  if(j>0 && j<3) continue;
501  fEFFproj[part][j] = nom->Projection(j,1,2,"E");
502  fEFFproj[part][j]->SetName(Form("NOMidx%dj%d12",part,j));
503 
504  TH1 *den=denom->Projection(j,1,2,"E");
505  den->SetName(Form("DENidx%dj%d12",part,j));
506  fEFFproj[part][j]->Divide( den );
507 
508  if(den) delete den;
509  }
510  break;
511  default:
512  break;
513  }
514 
515 
516  }
517 
518  // calculate main effiency map
519  nom->Divide(denom);
520  */
521 }
522 
523 Double_t CbmFastSim::GetEfficiency(Int_t idx, Double_t* vals, Int_t* coord) {
527 
528  if (fMethod == kIgnoreFluct) return 0.;
529 
530 
531  Int_t ndim = 0;
532  Double_t meanEff = (fMethod == kFactorize ? 1. : 0.);
533  Double_t eff = 0.;
534  // for (Int_t j = 0; j < 10; ++j) {
535  for (Int_t j = 10; j >= 0; j--) {
536 
537  if (!fEFFproj[idx][j]) continue;
538 
540  if (coord)
541  eff = fEFFproj[idx][j]->GetBinContent(coord[j]);
542  else
543  eff = fEFFproj[idx][j]->GetBinContent(fEFFproj[idx][j]->FindBin(vals[j]));
544  // Info("GetEfficiency","int eff for j:%d and idx:%d and value: %.2e is %.2e!",j,idx,vals[j],eff);
545 
547  switch (fMethod) {
548  case kIgnoreFluct: return 0.;
549  case kInterpolate: {
550  if (coord) {
553  for (Int_t c = 0; c < 3; c++) {
554  if (coord[c] == 0) coord[c] += 2;
555  if (coord[c] == 1) coord[c] += 1;
556  if (coord[c] == fEFF[idx]->GetAxis(c)->GetNbins()) coord[c] -= 2;
557  if (coord[c] == fEFF[idx]->GetAxis(c)->GetNbins() - 1)
558  coord[c] -= 1;
559  }
560  eff = fEFFproj[idx][j]->Interpolate(
561  fEFFproj[idx][j]->GetXaxis()->GetBinCenter(coord[0]),
562  fEFFproj[idx][j]->GetYaxis()->GetBinCenter(coord[1]),
563  fEFFproj[idx][j]->GetZaxis()->GetBinCenter(coord[2]));
564  } else
565  eff = fEFFproj[idx][j]->Interpolate(vals[0], vals[1], vals[2]);
566  return eff;
567  } break;
568  case kAverage:
569  case kLastDim:
570  if (eff < 1.e-10)
571  return 0.;
572  meanEff += eff;
573  ndim++;
574  break;
575  case kFactorize:
576  if (eff < 1.e-10)
577  return 0.;
578  meanEff *= eff;
579  ndim = 1; // ndim++;
580  break;
581  }
582  }
583 
584  // protect against division by 0
585  if (!ndim) {
586  Error("GetEfficiency", "number of dimensions for idx:%d is zero!", idx);
587  return 0.;
588  }
589 
590  // return mean efficiency
591  return (meanEff / ndim);
592 }
593 
594 Bool_t CbmFastSim::IsSelected(TObject* sel, Int_t opt) {
595  //
596  // check if partilce passed the cuts
597  //
598 
599  Int_t pdg = opt;
600  TVector3* mom = dynamic_cast<TVector3*>(sel);
601 
602  // Info("IsSelected","Compare pdg code %d to e(%d),mu(%d),pi(%d),K(%d),p(%d)",
603  // TMath::Abs(pdg),kElectron,kMuonMinus,kPiPlus,kKPlus,kProton);
604 
605  Int_t idx = -1;
606  switch (TMath::Abs(pdg)) {
607  case kGamma: idx = CbmFastSim::kGam; break;
608  case kElectron: idx = CbmFastSim::kEle; break;
609  case kMuonMinus: idx = CbmFastSim::kMuo; break;
610  case kPiPlus: idx = CbmFastSim::kPio; break;
611  case kKPlus: idx = CbmFastSim::kKao; break;
612  case kProton: idx = CbmFastSim::kPro; break;
613  default: /*Warning("IsSelected","Pdg code %d not supported, particle will not be processed",pdg);*/
614  return kFALSE;
615  }
616 
617  // Double_t rndm2 = fRand->Rndm();
618  // if(rndm2<0.1) return kTRUE;
619  // else return kFALSE;
620 
622  if (!fEFF[idx]) {
623  Warning("IsSelected", "No filter for pdg code %d", pdg);
624  return kTRUE;
625  }
626 
628  Double_t* vals = new Double_t[3];
629 
632  vals[0] = mom->Pt();
633  vals[1] = mom->Theta();
634  vals[2] = TVector2::Phi_0_2pi(mom->Phi());
635 
637  Long_t bin = fEFF[idx]->GetBin(vals, kFALSE);
638 
640  delete vals;
641 
642  Double_t eff = 0.0;
643  if (bin > 0) eff = fEFF[idx]->GetBinContent(bin);
644  Double_t rndm = fRand->Rndm();
645 
646 
647  if (rndm > eff)
648  return kFALSE;
649  else
650  return kTRUE;
651 }
652 
CbmFastSim::fRand
TRandom3 * fRand
Definition: CbmFastSim.h:65
CbmFastSim::Exec
virtual void Exec(Option_t *opt)
Definition: CbmFastSim.cxx:211
CbmFastSim::kPio
@ kPio
Definition: CbmFastSim.h:30
CbmFastSim::fMethod
EEfficiencyMethod fMethod
Definition: CbmFastSim.h:68
CbmFastSim::fgkNParts
static const Int_t fgkNParts
Definition: CbmFastSim.h:67
PairAnalysisHelper.h
CbmFastSim::kAverage
@ kAverage
Definition: CbmFastSim.h:35
CbmFastSim::fFastTracks
TClonesArray * fFastTracks
Definition: CbmFastSim.h:80
CbmFastSim::GetEfficiency
Double_t GetEfficiency(Int_t idx, Double_t *vals, Int_t *coord)
Definition: CbmFastSim.cxx:523
CbmFastSim::SetParContainers
virtual void SetParContainers()
Definition: CbmFastSim.cxx:184
CbmFastSim::kEle
@ kEle
Definition: CbmFastSim.h:30
CbmFastSim::kInterpolate
@ kInterpolate
Definition: CbmFastSim.h:34
CbmFastSim::fRFp
TH2F * fRFp[fgkNParts]
Definition: CbmFastSim.h:87
CbmFastSim::EParticleType
EParticleType
Definition: CbmFastSim.h:30
CbmFastSim::Smearing
Bool_t Smearing(TLorentzVector *p4, Int_t pdg)
Definition: CbmFastSim.cxx:365
CbmFastSim::Init
virtual InitStatus Init()
Definition: CbmFastSim.cxx:104
CbmFastSim::AnalyzeMap
Double_t AnalyzeMap(TH2F *map, Double_t refValue)
Definition: CbmFastSim.cxx:409
CbmFastSim::IsSelected
virtual Bool_t IsSelected(TObject *sel, Int_t opt)
Definition: CbmFastSim.cxx:594
CbmStack.h
CbmFastSim::Register
void Register()
Definition: CbmFastSim.cxx:90
CbmFastSim.h
CbmFastSim
Definition: CbmFastSim.h:27
CbmFastSim::kFactorize
@ kFactorize
Definition: CbmFastSim.h:36
CbmFastSim::PassEfficiencyFilter
Bool_t PassEfficiencyFilter(Int_t pdg, Double_t *vals, Int_t *coord)
Definition: CbmFastSim.cxx:305
CbmFastSim::~CbmFastSim
~CbmFastSim()
Definition: CbmFastSim.cxx:74
ClassImp
ClassImp(CbmConverterManager) InitStatus CbmConverterManager
Definition: CbmConverterManager.cxx:12
PairAnalysisHelper::CumulateSlicesX
void CumulateSlicesX(TH2 *h, Bool_t reverse=kFALSE, Bool_t norm=kFALSE)
Definition: PairAnalysisHelper.cxx:597
CbmStack::GetNtrack
virtual Int_t GetNtrack() const
Definition: CbmStack.h:179
CbmFastSim::fEFFdenom
THnBase * fEFFdenom[fgkNParts]
Definition: CbmFastSim.h:84
kGamma
@ kGamma
Definition: CbmLmvmHist.h:10
CbmFastSim::SetLookupEfficiency
void SetLookupEfficiency(THnBase *nom, THnBase *denom, EParticleType part)
Definition: CbmFastSim.cxx:455
CbmFastSim::fdbPdg
TDatabasePDG * fdbPdg
Definition: CbmFastSim.h:66
CbmFastSim::CbmFastSim
CbmFastSim(bool persist=true)
Definition: CbmFastSim.cxx:55
CbmFastSim::kLastDim
@ kLastDim
Definition: CbmFastSim.h:37
CbmFastSim::Finish
virtual void Finish()
Definition: CbmFastSim.cxx:203
CbmFastSim::kPro
@ kPro
Definition: CbmFastSim.h:30
CbmMCTrack.h
CbmFastSim::kMuo
@ kMuo
Definition: CbmFastSim.h:30
CbmStack::GetParticle
TParticle * GetParticle(Int_t trackId) const
Definition: CbmStack.cxx:420
CbmFastSim::InitTask
virtual void InitTask()
Definition: CbmFastSim.cxx:101
CbmFastSim::kIgnoreFluct
@ kIgnoreFluct
Definition: CbmFastSim.h:33
CbmFastSim::fEFF
THnBase * fEFF[fgkNParts]
Definition: CbmFastSim.h:83
CbmFastSim::fEFFproj
TH1 * fEFFproj[fgkNParts][10]
Definition: CbmFastSim.h:85
CbmStack
Definition: CbmStack.h:45
CbmFastSim::SetSeed
void SetSeed(unsigned int seed=65539)
Definition: CbmFastSim.cxx:195
CbmFastSim::kKao
@ kKao
Definition: CbmFastSim.h:30
CbmFastSim::kGam
@ kGam
Definition: CbmFastSim.h:30