PandaRoot
PndMdtParamDigi.h
Go to the documentation of this file.
1 //****************************************************************************
2 //* This file is part of PandaRoot. *
3 //* *
4 //* PandaRoot is distributed under the terms of the *
5 //* GNU General Public License (GPL) version 3, *
6 //* copied verbatim in the file "LICENSE". *
7 //* *
8 //* Copyright (C) 2006 - 2024 FAIR GmbH and copyright holders of PandaRoot *
9 //* The copyright holders are listed in the file "COPYRIGHTHOLDERS". *
10 //* The authors are listed in the file "AUTHORS". *
11 //****************************************************************************
12 
13 //==================================================================
14 // description
15 // parameterized digitizaiton model from garfield simulation
16 // author, Jifeng Hu, hu@to.infn.it, Torino University
17 //
18 
19 #ifndef PNDMDTPARAMDIGI_H
20 #define PNDMDTPARAMDIGI_H 1
21 
22 #include "TVector3.h"
23 #include "TGraphErrors.h"
24 #include <vector>
25 #include "TF1.h"
26 #include "TH1F.h"
27 #include "TH2F.h"
28 #include "TRandom.h"
29 #include <math.h>
30 #include <map>
31 
32 // value and its error
33 typedef std::pair<Double_t, Double_t> ValueErrorType;
34 
35 namespace CLHEP {
36 class RandGeneral;
37 }
38 
39 class PndMdtParamDigi : public TNamed {
40 
41  public:
44 
46  ~PndMdtParamDigi();
47  void SetVerbose(Int_t v) { fVerbose = v; }
48  // call Init before invoking the other interfaces.
49  Bool_t Init();
50  // signal on Anode wire
51  //==========================================================================
52  // Int_t particleType: electron, muon, pion, kaon+-, proton, indiced as 0, 1, 2, 3, 4
53  // TVector3 iniP : initial 3-momentum of this charged particle
54  // TVector3 iniPos : initial position when entering into a tube
55  // TVector3 fnlPos : final position when exiting a tube
56  //==========================================================================
57  // by default, strip lenght 100 cm
58  // this inteface should be invoked first when parameters change
59  PndMdtParamDigi &SetParams(Int_t ptlType, TVector3 iniP, TVector3 iniPos, TVector3 finalPos, Double_t stripLen = 100.);
60 
61  void UseNoise(Bool_t swith) { fUseNoise = swith; }
62  void UseDetailedSim(Bool_t swith = kTRUE) { fDetailedSim = swith; }
63  void UsePlot(Bool_t swith = kTRUE) { fUsePlot = swith; }
64  void UseGaussianAmp(Bool_t swith) { fGaussianAmp = swith; }
65  void SetOptimization(Int_t val) { fNumofTruncation = val; }
66  void SetNoiseWidth(Double_t anode, Double_t strip)
67  {
68  fNoiseSigmaAnode = anode;
69  fNoiseSigmaStrip = strip;
70  }
71  //===============================================================
72  // readout interfaces
73  // only fired information is provided
74  // for wire, <-1, fired time>, for strips, <strip number from 0, fired time>
75  std::vector<std::pair<Int_t, Double_t>> GetFiredInfo();
76  // both anode and strip signals
77  void Compute(Bool_t useConvolution = kTRUE);
78  // access data of induced current
79  const std::vector<Double_t> &GetWireSignal() const { return fSignalDataAnode; }
80  const std::map<Int_t, std::vector<Double_t>> &GetStripSignals() const { return fSignalDataStripM; }
81  //===============================================================
82 
83  // plot a charged track and its produced signal
84  void Draw();
85 
86  private:
87  // anode signal only
88  void GetSignal(Bool_t useConvolution = kTRUE);
89  // get threshold-crossing time and amplitude
90  // if tube is fired, return kTRUE;
91  // in case of kTRUE, time and amplitude will be set.
92  Bool_t Digitize(Double_t &time, Double_t &amp); // anode
93  Bool_t Digitize(Int_t stripNo, Double_t &time, Double_t &amp); // strip
94  //
95  void GetRawSignalbySimAvalanche(Double_t fNoiseLevel = 1.);
96  void GetRawSignalbyWeightingAvalanche(Double_t fNoiseLevel = 1.);
97  void AddNoise(Double_t fNoiseLevel = 1., Int_t isAnode = 1);
98  void ApplyTransferFunction(Double_t *fSignalData, Int_t nSize = fSamplingSize);
99 
100  // by default, sampling rate is set to 100 MHz, 200 sampling points
101  static const Int_t fSamplingSize = 20;
102  Int_t fSamplingRate;
103  Int_t fSamplingInterval;
104  std::vector<Double_t> fSignalDataAnode;
105  std::map<Int_t, std::vector<Double_t>> fSignalDataStripM; // strip no, and its corresponding signal
106  TF1 *fTrF; // transfer function
107  Bool_t fUseNoise;
108  Bool_t fGaussianAmp; // amplication model
109  Double_t fNoiseSigmaAnode; // 0.01 mu Ampere
110  Double_t fNoiseSigmaStrip; // 0.01 mu Ampere
111  Int_t fNumofTruncation;
112 
113  private:
114  // drift time of primary ionized electrons to wire surface
115  inline Double_t GetElectronDriftTime(TVector2 iniPos)
116  {
117  const Double_t p0 = 5.35513e+03;
118  const Double_t p1 = 3.92749e+02;
119  const Double_t p2 = -5.41433e+03;
120  const Double_t p3 = -3.31810e+03;
121  Double_t x = iniPos.Mod() - cWIRERADIUS;
122  Double_t x1 = log(1. + x);
123  return (x * (p0 + p1 * pow(x, 0.25)) + x1 * (p2 + p3 * x1));
124  }
125  // drift time of secondary ions from wire surface to a given point
126  inline Double_t GetShiftTime(TVector2 fWireSurfacePos, TVector2 fProductionPos)
127  {
128  const Double_t p0 = 1.81940e-8; // unit: micro sencond
129  const Double_t p1 = 1.80248;
130  const Double_t p2 = 3.68652e+2;
131  const Double_t p3 = 1.67353e+3;
132  // const Double_t mean = 0; //[R.K. 01/2017] unused variable?
133  // const Double_t sigma = 2.73477e-3;//micro second //[R.K. 01/2017] unused variable?
134  Double_t dr = (fProductionPos - fWireSurfacePos).Mod();
135  return (((p3 * dr + p2) * dr + p1) * dr + p0 /*+ gRandom->Gaus(mean, sigma)*/) * 1.e3; // nano sencod
136  }
137 
138  private:
139  Int_t fParticleType; // electron 0, muon 1, pion 2,
140  TVector3 fParticleMomentum; // momentum of the incident charged particle
141  TVector3 fInitPostion;
142  TVector3 fExitPosition;
143  Bool_t fParamsChanged;
144  Double_t fStripLength;
145  Double_t fRestMass[5];
146 
147  void GetMPVofPrimaryIonization(Int_t particleType, const TVector3 &momentum, ValueErrorType &val) const;
148  Int_t GetAmplicationFactor(Int_t particleType, Double_t momentum) const;
149  // sampling position of secondary ion/electron for a given center
150  Bool_t fDetailedSim;
151  void SamplingPosition(TVector2 fDirection, TVector2 &fIonProductionPos);
152  // mean free path function as velocity
153  TGraphErrors *gFreePath[5]; //
154 
155  inline Double_t GetMeanFreePath(Int_t ptlType, Double_t mom) const
156  {
157  Double_t mMass = fRestMass[ptlType];
158  // Double_t mEnergy = sqrt(mom*mom + mMass*mMass);
159  Double_t mBeta = mom / mMass;
160  Double_t mPath = gFreePath[ptlType]->Eval(mBeta);
161  return gRandom->Exp(mPath);
162  }
163 
164  // avalanche size function
165  Bool_t fUseAvaHist;
166  TH1F *hAvaSize;
167  // TH1F* hAnodeI;
168  // TH1F* hCathodeIx;
169  // TH1F* hCathodeIy;
170  //
172  static const Int_t NRAD = 100; // radial bins, not changeable
173  static const Int_t NPHI = 100; // phi bins, not changeable
174  // geometry constant
175  Double_t cWIRERADIUS; // unit: cm
176  Double_t cCELLSIZE; // unit: cm
177  Double_t cSTRIPWIDTH; // unit: cm
178  Double_t cMAXRADIUS; // not changeable
179  Double_t cMAPFACTOR; // not changeable
180  Double_t cUNITLOGDR; //= TMath::Log(1. + cMAPFACTOR*(cMAXRADIUS - cWIRERADIUS))/NRAD;
181  Double_t cUNITPHI; //= TMath::Pi()*2./NPHI;
182  TH1F *hAnodeI1d[NRAD];
183  TH1F *hCathodeI2d[NRAD][NPHI];
184  Int_t fNr[NRAD];
185  Int_t fNrNphi[NRAD * NPHI];
186 
187  // cluster initialization
188  struct ClusInfo {
189  explicit ClusInfo(TVector3 _v = TVector3()) : fPosition(_v) {}
190  TVector3 fPosition;
191  Int_t fNumofPrimaryIonization;
192  Int_t fNumofIonsofThisCluster;
193  Double_t fStartTime;
194  };
195 
196  // histograms
197  Bool_t fUsePlot;
198  TH1F *hIndex;
199  TH2F *h2dIndex;
200  TH2F *hTrack;
201  TH1F *hIonNum;
202  TH2F *hAvaPos;
203  TH1F *hAvaNum;
204  TH1F *hAva1D;
205  TH2F *hAva2D;
206  //
207  TH1F *hExp1;
208  TH1F *hExp2;
209  TH1F *hGen;
210  TH1F *hLan;
211  Int_t fVerbose;
212  //
213  public:
214  struct AvaBinType {
215  AvaBinType(Int_t _i = 0, Double_t _v = 0.) : Index(_i), Probabilty(_v) {}
216  Int_t Index;
217  Double_t Probabilty;
218  }; // index, density
219  private:
220  std::vector<AvaBinType> fProbFunc1D;
221  std::vector<AvaBinType> fProbFunc2D;
222 
223  //
224  Double_t SPEEDOFLIGHT;
225  CLHEP::RandGeneral *fRandAva;
226 
227  ClassDef(PndMdtParamDigi, 1);
228 };
229 
230 #endif
std::pair< Double_t, Double_t > ValueErrorType
void UseGaussianAmp(Bool_t swith)
void SetNoiseWidth(Double_t anode, Double_t strip)
friend F32vec4 log(const F32vec4 &a)
Definition: P4_F32vec4.h:129
__m128 v
Definition: P4_F32vec4.h:15
const std::vector< Double_t > & GetWireSignal() const
void UsePlot(Bool_t swith=kTRUE)
const std::map< Int_t, std::vector< Double_t > > & GetStripSignals() const
AvaBinType(Int_t _i=0, Double_t _v=0.)
void UseDetailedSim(Bool_t swith=kTRUE)
void UseNoise(Bool_t swith)
void SetOptimization(Int_t val)
void SetVerbose(Int_t v)