BOSS 7.1.2
BESIII Offline Software System
Loading...
Searching...
No Matches
BesTofDigitizerBrV1.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: BesTofDigitizerBrV1.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
26
30
32{
33 G4cout<<"BesTofDigitizerBrV1::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 partId = hit->GetPartId();
78 G4int scinNb = hit->GetScinNb();
79 G4double time = hit->GetTime();
81 {
83 m_trackIndex = trackIndex;
84 }
85 G4double edep = hit->GetEdep();
86 G4ThreeVector pos = hit->GetPos();
87 G4double posx = pos.x();
88 G4double posy = pos.y();
89 G4double posz = pos.z();
90
91 G4double atten;
92 atten = m_tofCaliSvc->BAtten(scinNb);
93 G4double pathL[2];
94 pathL[0]=m_scinLength/2-posz;
95 pathL[1]=posz+m_scinLength/2;
96 for (G4int j=0;j<2;j++)
97 {
98 //G4cout<<"atten:"<<atten<<G4endl;
99 m_ADC[j] += edep*exp(-pathL[j]/atten);
100 }
101 if (time<m_t1st)
102 {
103 m_t1st = time;
104 m_z = posz;
105 m_TDC[0] = m_TDC[1] = m_t1st;
106 }
107}
108
110{
111 G4cout<<"m_t1st:"<<m_t1st<<" m_z:"<<m_z<<G4endl;
112
113 /*G4double tofRes = 0.08;
114 G4double atten;
115 for(G4int i=0;i<2;i++)
116 {
117 m_TDC[i] += tofRes * G4RandGauss::shoot();
118 }*/
119
120 double pp1[10];
121 double pp2[10];
122 for (int i=0;i<10;i++)
123 {
124 pp1[i]=m_tofCaliSvc->BTof(scinNb)->getP1(i);
125 pp2[i]=m_tofCaliSvc->BTof(scinNb)->getP2(i);
126 //G4cout<<"pp1["<<i<<"]:"<<pp1[i]<<" "<<"pp2["<<i<<"]:"<<pp2[i]<<G4endl;
127 }
128 double pp[10];
129 for (int i=0;i<2;i++)
130 {
131 m_ADC[i] *= 7.;
132 for (int j=0;j<10;j++)
133 {
134 if (i==0)
135 pp[j]=pp1[j];
136 else
137 pp[j]=pp2[j];
138 }
139 m_TDC[i] += (pp[0]+pp[1]*m_z)/TMath::Sqrt(m_ADC[i])+
140 pp[2]*m_ADC[i]+
141 pp[3]*m_ADC[i]*m_ADC[i]+
142 pp[4]*m_ADC[i]*m_ADC[i]*m_ADC[i]+
143 pp[5]/(84.2*84.2+m_z*m_z)+
144 pp[6]*m_z+
145 pp[7]*m_z*m_z+
146 pp[8]*m_z*m_z*m_z+
147 pp[9];
148 }
149}
G4TDigiCollection< BesTofDigi > BesTofDigitsCollection
Definition BesTofDigi.hh:79
G4THitsCollection< BesTofHit > BesTofHitsCollection
Definition BesTofHit.hh:116
Double_t time
EvtComplex exp(const EvtComplex &c)
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
void TofPmtAccum(BesTofHit *)
virtual void Digitize(ScintSingle *, BesTofDigitsCollection *)
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 GetPartId()
Definition BesTofHit.hh:68
G4int GetTrackIndex()
Definition BesTofHit.hh:66
virtual BTofCal * BTof(unsigned id) const =0
virtual const double BAtten(unsigned id)=0
G4int GetPartId()
vector< G4int > * GetHitIndexes()
G4int GetScinNb()