BOSS 8.0.0
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
12#include "TofSim/BesTofDigitizerBrV1.hh"
13#include "G4DigiManager.hh"
14#include "Randomize.hh"
15#include "TofSim/BesTofDigi.hh"
16#include "TofSim/BesTofGeoParameter.hh"
17#include "TofSim/BesTofHit.hh"
18#include <TMath.h>
19#include <math.h>
20
25
27
29 G4cout << "BesTofDigitizerBrV1::Digitize" << G4endl;
31 G4DigiManager* digiManager = G4DigiManager::GetDMpointer();
32 G4int THCID = digiManager->GetHitsCollectionID( "BesTofHitsCollection" );
33 m_THC = (BesTofHitsCollection*)( digiManager->GetHitsCollection( THCID ) );
34 if ( m_THC )
35 {
36 G4int partId, scinNb, nHits;
37 BesTofHit* hit;
38 partId = scint->GetPartId();
39 scinNb = scint->GetScinNb();
40 nHits = scint->GetHitIndexes()->size();
41 TofPmtInit();
42 for ( G4int j = 0; j < nHits; j++ )
43 {
44 hit = ( *m_THC )[( *( scint->GetHitIndexes() ) )[j]];
45 TofPmtAccum( hit );
46 }
47
48 Smear( scinNb );
49 if ( m_TDC[0] > 0 )
50 {
51 BesTofDigi* digi = new BesTofDigi;
53 digi->SetPartId( partId );
54 digi->SetScinNb( scinNb );
55 digi->SetForwADC( m_ADC[0] );
56 digi->SetForwTDC( m_TDC[0] );
57 digi->SetBackADC( m_ADC[1] );
58 digi->SetBackTDC( m_TDC[1] );
59 m_besTofDigitsCollection->insert( digi );
60 }
61 }
62}
63
65 Initialize();
66 m_t1st = 9999.;
67}
68
70 G4int trackIndex = hit->GetTrackIndex();
71 G4int partId = hit->GetPartId();
72 G4int scinNb = hit->GetScinNb();
73 G4double time = hit->GetTime();
74 if ( time < m_globalTime )
75 {
77 m_trackIndex = trackIndex;
78 }
79 G4double edep = hit->GetEdep();
80 G4ThreeVector pos = hit->GetPos();
81 G4double posx = pos.x();
82 G4double posy = pos.y();
83 G4double posz = pos.z();
84
85 G4double atten;
86 atten = m_tofCaliSvc->BAtten( scinNb );
87 G4double pathL[2];
88 pathL[0] = m_scinLength / 2 - posz;
89 pathL[1] = posz + m_scinLength / 2;
90 for ( G4int j = 0; j < 2; j++ )
91 {
92 // G4cout<<"atten:"<<atten<<G4endl;
93 m_ADC[j] += edep * exp( -pathL[j] / atten );
94 }
95 if ( time < m_t1st )
96 {
97 m_t1st = time;
98 m_z = posz;
99 m_TDC[0] = m_TDC[1] = m_t1st;
100 }
101}
102
103void BesTofDigitizerBrV1::Smear( G4int scinNb ) {
104 G4cout << "m_t1st:" << m_t1st << " m_z:" << m_z << G4endl;
105
106 /*G4double tofRes = 0.08;
107 G4double atten;
108 for(G4int i=0;i<2;i++)
109 {
110 m_TDC[i] += tofRes * G4RandGauss::shoot();
111 }*/
112
113 double pp1[10];
114 double pp2[10];
115 for ( int i = 0; i < 10; i++ )
116 {
117 pp1[i] = m_tofCaliSvc->BTof( scinNb )->getP1( i );
118 pp2[i] = m_tofCaliSvc->BTof( scinNb )->getP2( i );
119 // G4cout<<"pp1["<<i<<"]:"<<pp1[i]<<" "<<"pp2["<<i<<"]:"<<pp2[i]<<G4endl;
120 }
121 double pp[10];
122 for ( int i = 0; i < 2; i++ )
123 {
124 m_ADC[i] *= 7.;
125 for ( int j = 0; j < 10; j++ )
126 {
127 if ( i == 0 ) pp[j] = pp1[j];
128 else pp[j] = pp2[j];
129 }
130 m_TDC[i] += ( pp[0] + pp[1] * m_z ) / TMath::Sqrt( m_ADC[i] ) + pp[2] * m_ADC[i] +
131 pp[3] * m_ADC[i] * m_ADC[i] + pp[4] * m_ADC[i] * m_ADC[i] * m_ADC[i] +
132 pp[5] / ( 84.2 * 84.2 + m_z * m_z ) + pp[6] * m_z + pp[7] * m_z * m_z +
133 pp[8] * m_z * m_z * m_z + pp[9];
134 }
135}
Double_t time
EvtComplex exp(const EvtComplex &c)
G4TDigiCollection< BesTofDigi > BesTofDigitsCollection
G4THitsCollection< BesTofHit > BesTofHitsCollection
void TofPmtAccum(BesTofHit *)
virtual void Digitize(ScintSingle *, BesTofDigitsCollection *)
static BesTofGeoParameter * GetInstance()