BOSS 7.1.0
BESIII Offline Software System
Loading...
Searching...
No Matches
BesTofDigitizerEcV1.cc
Go to the documentation of this file.
1//---------------------------------------------------------------------------//
2// BOOST --- BESIII Object_Oriented Simulation Tool //
3//---------------------------------------------------------------------------//
4//Description:
5//Author: Dengzy
6//Created: Mar, 2004
7//Modified:
8//Comment:
9//---------------------------------------------------------------------------//
10//$Id: BesTofDigitizerEcV1.cc
11
13#include "BesTofDigi.hh"
14#include "BesTofHit.hh"
15#include "G4DigiManager.hh"
16#include "BesTofGeoParameter.hh"
17#include "Randomize.hh"
18#include "TMath.h"
19#include <math.h>
20
22{
24 m_bucketPosR = tofPara->GetBucketPosR(); // 445 ???
25}
26
28{
29}
30
32{
33 G4cout<<"BesTofDigitizerEcV1::Digitize"<<G4endl;
35 G4DigiManager* digiManager = G4DigiManager::GetDMpointer();
36 G4int THCID = digiManager->GetHitsCollectionID("BesTofHitsCollection");
37 m_THC = (BesTofHitsCollection*) (digiManager->GetHitsCollection(THCID));
38 if (m_THC)
39 {
40 G4int partId, scinNb, nHits;
41 BesTofHit* hit;
42 partId=scint->GetPartId();
43 scinNb=scint->GetScinNb();
44 nHits=scint->GetHitIndexes()->size();
45 TofPmtInit();
46 for (G4int j=0;j<nHits;j++)
47 {
48 hit= (*m_THC)[( *(scint->GetHitIndexes()) )[j]];
49 TofPmtAccum(hit);
50 }
51
52 Smear(scinNb);
53 if ( m_TDC[0]>0 )
54 {
55 BesTofDigi* digi = new BesTofDigi;
57 digi->SetPartId(partId);
58 digi->SetScinNb(scinNb);
59 digi->SetForwADC( m_ADC[0]) ;
60 digi->SetForwTDC( m_TDC[0]) ;
61 digi->SetBackADC( m_ADC[1]) ;
62 digi->SetBackTDC( m_TDC[1]) ;
63 m_besTofDigitsCollection->insert(digi);
64 }
65 }
66}
67
69{
70 Initialize();
71 m_t1st = 9999.;
72}
73
75{
76 G4int trackIndex = hit->GetTrackIndex();
77 G4int scinNb = hit->GetScinNb();
78 G4double time = hit->GetTime();
80 {
82 m_trackIndex = trackIndex;
83 }
84 G4double edep = hit->GetEdep();
85 G4ThreeVector pos = hit->GetPos();
86 G4double posx = pos.x();
87 G4double posy = pos.y();
88 G4double pathL=abs(m_bucketPosR-sqrt(posx*posx+posy*posy));
89 G4double atten;
90 atten = m_tofCaliSvc->EAtten(scinNb);
91 m_ADC[0] += edep*exp(-pathL/atten);
92
93 if (time<m_t1st)
94 {
95 m_t1st = time;
96 m_r = sqrt(posx*posx+posy*posy);
97 m_TDC[0] = m_t1st;
98 }
99}
100
102{
103 /*G4double tofRes = 0.08;
104 for(G4int i=0;i<2;i++)
105 {
106 m_TDC[i] += tofRes * G4RandGauss::shoot();
107 }*/
108
109 double pp[8];
110 for (int i=0;i<8;i++)
111 {
112 pp[i]=m_tofCaliSvc->ETof(scinNb)->getP(i);
113 }
114 m_ADC[0] *= 7.;
115 m_TDC[0] += (pp[0]+pp[1]*m_r)/TMath::Sqrt(m_ADC[0])+
116 pp[2]/m_ADC[0]+
117 pp[3]*m_r/m_ADC[0]+
118 pp[4]*m_r+
119 pp[5]*m_r*m_r+
120 pp[6]*m_r*m_r*m_r+
121 pp[7];
122}
G4TDigiCollection< BesTofDigi > BesTofDigitsCollection
Definition: BesTofDigi.hh:79
G4THitsCollection< BesTofHit > BesTofHitsCollection
Definition: BesTofHit.hh:116
Double_t time
EvtComplex exp(const EvtComplex &c)
Definition: EvtComplex.hh:252
void SetPartId(G4int partId)
Definition: BesTofDigi.hh:37
void SetForwADC(G4double ADC)
Definition: BesTofDigi.hh:40
void SetTrackIndex(G4int index)
Definition: BesTofDigi.hh:36
void SetBackADC(G4double ADC)
Definition: BesTofDigi.hh:41
void SetForwTDC(G4double TDC)
Definition: BesTofDigi.hh:42
void SetScinNb(G4int scinNb)
Definition: BesTofDigi.hh:39
void SetBackTDC(G4double TDC)
Definition: BesTofDigi.hh:43
virtual void Digitize(ScintSingle *, BesTofDigitsCollection *)
void TofPmtAccum(BesTofHit *)
ITofCaliSvc * m_tofCaliSvc
BesTofDigitsCollection * m_besTofDigitsCollection
BesTofHitsCollection * m_THC
static BesTofGeoParameter * GetInstance()
G4double GetTime()
Definition: BesTofHit.hh:74
G4double GetEdep()
Definition: BesTofHit.hh:70
G4ThreeVector GetPos()
Definition: BesTofHit.hh:73
G4int GetScinNb()
Definition: BesTofHit.hh:69
G4int GetTrackIndex()
Definition: BesTofHit.hh:66
virtual ETofCal * ETof(unsigned id) const =0
virtual const double EAtten(unsigned id)=0
G4int GetPartId()
Definition: ScintSingle.hh:44
vector< G4int > * GetHitIndexes()
Definition: ScintSingle.hh:47
G4int GetScinNb()
Definition: ScintSingle.hh:45