BOSS 8.0.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
12#include "TofSim/BesTofDigitizerEcV1.hh"
13#include "G4DigiManager.hh"
14#include "Randomize.hh"
15#include "TMath.h"
16#include "TofSim/BesTofDigi.hh"
17#include "TofSim/BesTofGeoParameter.hh"
18#include "TofSim/BesTofHit.hh"
19#include <math.h>
20
23 m_bucketPosR = tofPara->GetBucketPosR(); // 445 ???
24}
25
27
29 G4cout << "BesTofDigitizerEcV1::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 scinNb = hit->GetScinNb();
72 G4double time = hit->GetTime();
73 if ( time < m_globalTime )
74 {
76 m_trackIndex = trackIndex;
77 }
78 G4double edep = hit->GetEdep();
79 G4ThreeVector pos = hit->GetPos();
80 G4double posx = pos.x();
81 G4double posy = pos.y();
82 G4double pathL = abs( m_bucketPosR - sqrt( posx * posx + posy * posy ) );
83 G4double atten;
84 atten = m_tofCaliSvc->EAtten( scinNb );
85 m_ADC[0] += edep * exp( -pathL / atten );
86
87 if ( time < m_t1st )
88 {
89 m_t1st = time;
90 m_r = sqrt( posx * posx + posy * posy );
91 m_TDC[0] = m_t1st;
92 }
93}
94
95void BesTofDigitizerEcV1::Smear( G4int scinNb ) {
96 /*G4double tofRes = 0.08;
97 for(G4int i=0;i<2;i++)
98 {
99 m_TDC[i] += tofRes * G4RandGauss::shoot();
100 }*/
101
102 double pp[8];
103 for ( int i = 0; i < 8; i++ ) { pp[i] = m_tofCaliSvc->ETof( scinNb )->getP( i ); }
104 m_ADC[0] *= 7.;
105 m_TDC[0] += ( pp[0] + pp[1] * m_r ) / TMath::Sqrt( m_ADC[0] ) + pp[2] / m_ADC[0] +
106 pp[3] * m_r / m_ADC[0] + pp[4] * m_r + pp[5] * m_r * m_r +
107 pp[6] * m_r * m_r * m_r + pp[7];
108}
Double_t time
EvtComplex exp(const EvtComplex &c)
G4TDigiCollection< BesTofDigi > BesTofDigitsCollection
G4THitsCollection< BesTofHit > BesTofHitsCollection
virtual void Digitize(ScintSingle *, BesTofDigitsCollection *)
void TofPmtAccum(BesTofHit *)
static BesTofGeoParameter * GetInstance()