BOSS 7.1.0
BESIII Offline Software System
Loading...
Searching...
No Matches
BesTofDigitizerBrV2 Class Reference

#include <BesTofDigitizerBrV2.hh>

+ Inheritance diagram for BesTofDigitizerBrV2:

Public Member Functions

 BesTofDigitizerBrV2 ()
 
 ~BesTofDigitizerBrV2 ()
 
virtual void Digitize (ScintSingle *, BesTofDigitsCollection *)
 
void ReadData ()
 
void TofPmtInit ()
 
void TofPmtAccum (BesTofHit *, G4int)
 
void DirectPh (G4int, G4ThreeVector, G4double &, G4int &)
 
G4double Scintillation (G4int)
 
G4double TransitTime ()
 
void AccuSignal (G4double, G4int)
 
void TofPmtRspns (G4int)
 
G4double Reflectivity (G4double n1, G4double n2, G4double theta)
 
G4double Reflectivity (G4double n1, G4double n2, G4double n3, G4double theta)
 
G4double BirksLaw (BesTofHit *hit)
 
- Public Member Functions inherited from BesTofDigitizerV
 BesTofDigitizerV ()
 
 ~BesTofDigitizerV ()
 
void Initialize ()
 
virtual void Digitize (ScintSingle *, BesTofDigitsCollection *)
 

Additional Inherited Members

- Protected Attributes inherited from BesTofDigitizerV
BesTofDigitsCollectionm_besTofDigitsCollection
 
BesTofHitsCollectionm_THC
 
ITofCaliSvcm_tofCaliSvc
 
ITofSimSvcm_tofSimSvc
 
ITofQElecSvcm_tofQElecSvc
 
G4double m_ADC [2]
 
G4double m_TDC [2]
 
G4int m_trackIndex
 
G4double m_globalTime
 
- Static Protected Attributes inherited from BesTofDigitizerV
static bool m_booked = false
 
static NTuple::Tuple * m_tupleTof1 = 0
 
static NTuple::Item< double > m_partId
 
static NTuple::Item< double > m_scinNb
 
static NTuple::Item< double > m_edep
 
static NTuple::Item< double > m_nHits
 
static NTuple::Item< double > m_time1st0
 
static NTuple::Item< double > m_time1st1
 
static NTuple::Item< double > m_timelast0
 
static NTuple::Item< double > m_timelast1
 
static NTuple::Item< double > m_totalPhot0
 
static NTuple::Item< double > m_totalPhot1
 
static NTuple::Item< double > m_NphAllSteps
 
static NTuple::Item< double > m_max0
 
static NTuple::Item< double > m_max1
 
static NTuple::Item< double > m_tdc0
 
static NTuple::Item< double > m_adc0
 
static NTuple::Item< double > m_tdc1
 
static NTuple::Item< double > m_adc1
 
static NTuple::Tuple * m_tupleTof2 = 0
 
static NTuple::Item< double > m_eTotal
 
static NTuple::Item< double > m_nDigi
 
static NTuple::Item< double > m_partIdMPV
 
static NTuple::Item< double > m_scinNbMPV
 
static NTuple::Item< double > m_edepMPV
 
static NTuple::Item< double > m_nDigiOut
 
static NTuple::Tuple * m_tupleTof3 = 0
 
static NTuple::Item< double > m_forb
 
static NTuple::Item< double > m_timeFlight
 
static NTuple::Item< double > m_ddT
 
static NTuple::Item< double > m_scinSwim
 
static NTuple::Item< double > m_scinTime
 
static NTuple::Item< double > m_transitTime
 
static NTuple::Item< double > m_endTime
 
static NTuple::Item< double > m_edepHit
 

Detailed Description

Definition at line 34 of file BesTofDigitizerBrV2.hh.

Constructor & Destructor Documentation

◆ BesTofDigitizerBrV2()

BesTofDigitizerBrV2::BesTofDigitizerBrV2 ( )

Definition at line 29 of file BesTofDigitizerBrV2.cc.

30{
31 ReadData();
32 m_timeBinSize=0.005;
33
34 //retrieve G4Svc
35 ISvcLocator* svcLocator = Gaudi::svcLocator();
36 IG4Svc* tmpSvc;
37 StatusCode sc = svcLocator->service("G4Svc", tmpSvc);
38 if(!sc.isSuccess())
39 {
40 std::cout << " Could not initialize Realization Service in BesTofDigitizerBrV2" << std::endl;
41 }
42 else
43 {
44 m_G4Svc = dynamic_cast<G4Svc *>(tmpSvc);
45 }
46
47 //retrieve RealizationSvc
48 IRealizationSvc *tmpReal;
49 StatusCode scReal = svcLocator->service("RealizationSvc",tmpReal);
50 if (!scReal.isSuccess())
51 {
52 std::cout << " Could not initialize Realization Service in BesTofDigitizerBrV2" << std::endl;
53 }
54 else
55 {
56 m_RealizationSvc = dynamic_cast<RealizationSvc*>(tmpReal);
57 }
58
59
60}
Definition: G4Svc.h:33
Definition: IG4Svc.h:31

◆ ~BesTofDigitizerBrV2()

BesTofDigitizerBrV2::~BesTofDigitizerBrV2 ( )

Definition at line 90 of file BesTofDigitizerBrV2.cc.

91{
92 ;
93}

Member Function Documentation

◆ AccuSignal()

void BesTofDigitizerBrV2::AccuSignal ( G4double  endTime,
G4int  forb 
)

Definition at line 752 of file BesTofDigitizerBrV2.cc.

753{
754 G4int ihst;
755 ihst=G4int(endTime/m_timeBinSize);
756 if (ihst>0 &&ihst<m_profBinN)
757 {
758 m_nPhot[ihst][forb]=m_nPhot[ihst][forb]+1;
759 m_totalPhot[forb]=m_totalPhot[forb]+1;
760 }
761}
const G4int m_profBinN

Referenced by TofPmtAccum().

◆ BirksLaw()

G4double BesTofDigitizerBrV2::BirksLaw ( BesTofHit hit)

Definition at line 385 of file BesTofDigitizerBrV2.cc.

386{
387 const G4double kappa = 0.015*cm/MeV;
388 const G4String brMaterial = "BC408";
389 G4double dE = hit->GetEdep();
390 //G4cout << "The edep is "<< dE << G4endl;
391 G4double dX = hit->GetStepL();
392 //G4Material* materiral = hit->GetMaterial();
393 G4double charge = hit->GetCharge();
394 G4double cor_dE = dE;
395 //if((materiral->GetName()==brMaterial) && charge!=0.&& dX!=0.)
396 if(charge!=0.&& dX!=0.)
397 {
398 cor_dE = dE/(1+kappa*dE/dX);
399 //if(dE>20)
400 //{
401 // G4cout << "\n dE > 20. Details are below:" << G4endl;
402 // G4cout << "dE/dx:" << dE/dX << G4endl;
403 // G4cout << "dE:" << dE << "; dX:" << dX << G4endl;
404 // G4cout << "It is BC408. cor_dE is " << cor_dE << G4endl;
405 // G4double ratio = cor_dE/dE;
406 // G4cout << "The ratio cor_dE/edep is "<< ratio << G4endl;
407 //}
408 //G4cout << "It is BC408. cor_dE is " << cor_dE << G4endl;
409 //G4double ratio = cor_dE/dE;
410 //G4cout << "The ratio cor_dE/edep is "<< ratio << G4endl;
411 }
412 return cor_dE;
413
414}
G4double GetStepL()
Definition: BesTofHit.hh:71
G4double GetCharge()
Definition: BesTofHit.hh:80
G4double GetEdep()
Definition: BesTofHit.hh:70
float charge

Referenced by TofPmtAccum().

◆ Digitize()

void BesTofDigitizerBrV2::Digitize ( ScintSingle scint,
BesTofDigitsCollection DC 
)
virtual

Reimplemented from BesTofDigitizerV.

Definition at line 95 of file BesTofDigitizerBrV2.cc.

96{
97 m_beamTime = m_G4Svc->GetBeamTime() * ns;
99
100 G4DigiManager* digiManager = G4DigiManager::GetDMpointer();
101 G4int THCID = digiManager->GetHitsCollectionID("BesTofHitsCollection");
102 m_THC = (BesTofHitsCollection*) (digiManager->GetHitsCollection(THCID));
103
104 if (m_G4Svc->TofRootFlag())
105 {
106 m_eTotal = 0;
107 m_nDigi = 0;
108 m_partIdMPV = -9;
109 m_scinNbMPV = -9;
110 m_edepMPV = 0;
111 m_nDigiOut = 0;
112 }
113
114 if (m_THC)
115 {
116 //for each digi, compute TDC and ADC
117 G4int partId, scinNb, nHits;
118 G4double edep;
119 BesTofHit* hit;
120 partId=scint->GetPartId();
121 scinNb=scint->GetScinNb();
122 edep = scint->GetEdep();
123 nHits=scint->GetHitIndexes()->size();
124
125
126 //std::cout << "BesTofDigitizerBrV2 Partid scinNb " << partId << " " << scinNb << std::endl;
127 //cout << "*** scinNb:" << scinNb << " *** " << m_tofCaliSvc->BAtten(scinNb) << "***" << endl;
128 //cout << "*****scinNb:"<< scinNb << " ***** A2:" << m_tofCaliSvc->BGainBackward(scinNb)
129 // << " ***** A1:" << m_tofCaliSvc->BGainForward(scinNb) << endl;
130
131 TofPmtInit();
132
133 //fill tof Ntuple
134 if (m_G4Svc->TofRootFlag())
135 {
136 if (edep>m_edepMPV)
137 {
138 m_partIdMPV = partId;
139 m_scinNbMPV = scinNb;
140 m_edepMPV = edep;
141 }
142 m_eTotal += edep;
143 m_nDigi ++;
144
145 m_partId = partId;
146 m_scinNb = scinNb;
147 m_edep = edep;
148 m_nHits = nHits;
149 }
150
151 if (edep>0.01)
152 {
153 for (G4int j=0;j<nHits;j++)
154 {
155 hit= (*m_THC)[( *(scint->GetHitIndexes()) )[j]];
156 TofPmtAccum(hit, scinNb);
157 }
158 if (m_G4Svc->TofRootFlag())
159 {
160 m_time1st0=m_t1st[0];
161 m_time1st1=m_t1st[1];
162 m_timelast0=m_tLast[0];
163 m_timelast1=m_tLast[1];
164 m_totalPhot0=m_totalPhot[0];
165 m_totalPhot1=m_totalPhot[1];
166 }
167 //get final tdc and adc
168 TofPmtRspns(scinNb);
169 //G4cout<<"pre-cut " << partId << "\nadc0:"<<m_ADC[0]<<"; adc1:"
170 // <<m_ADC[1]<<"; tdc0:"<<m_TDC[0]<<"; tdc1:"
171 // <<m_TDC[1]<<G4endl;
172
173 G4double temp0 = m_ADC[0]+m_TDC[0];
174 G4double temp1 = m_ADC[1]+m_TDC[1];
175 //cout << "partid: " << partId << " temp0: " << temp0 << " temp1: " << temp1 << endl;
176 //if ( partId==1 && m_ADC[0]>255 && m_ADC[1]>255 && m_TDC[0]>0. && m_TDC[1]>0.)
177 if ( partId==1 && (temp0 > -1998. || temp1 > -1998.))
178 {
179 //const double MAX_ADC = 8191 * 0.3; // channel set up to 8192 will lead to overflow.
180 BesTofDigi* digi = new BesTofDigi;
182 digi->SetPartId(partId);
183 digi->SetScinNb(scinNb);
184 digi->SetForwADC( m_ADC[0]) ;
185 digi->SetBackADC( m_ADC[1]) ;
186 if (m_TDC[0]>0)
187 m_TDC[0] = m_TDC[0]+m_beamTime;
188 if (m_TDC[1]>0)
189 m_TDC[1] = m_TDC[1]+m_beamTime;
190 digi->SetForwTDC( m_TDC[0]) ;
191 digi->SetBackTDC( m_TDC[1]) ;
192 //G4cout<<"+++++++++++++++++++++++++++++barrel\nadc0:"<<m_ADC[0]<<"; adc1:"
193 // <<m_ADC[1]<<"; tdc0:"<<m_TDC[0]<<"; tdc1:"
194 // <<m_TDC[1]<<G4endl;
195
196 m_besTofDigitsCollection->insert(digi);
197
198 if (m_G4Svc->TofRootFlag())
199 m_nDigiOut++;
200
201 }
202 if (m_G4Svc->TofRootFlag())
203 m_tupleTof1->write();
204
205 }
206
207 if (m_G4Svc->TofRootFlag())
208 m_tupleTof2->write();
209 }
210}
G4THitsCollection< BesTofHit > BesTofHitsCollection
Definition: BesTofHit.hh:116
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 *, G4int)
static NTuple::Item< double > m_partId
static NTuple::Item< double > m_eTotal
static NTuple::Item< double > m_timelast1
static NTuple::Item< double > m_timelast0
static NTuple::Item< double > m_scinNb
static NTuple::Tuple * m_tupleTof2
static NTuple::Item< double > m_edep
static NTuple::Item< double > m_totalPhot0
static NTuple::Item< double > m_scinNbMPV
static NTuple::Item< double > m_nDigi
BesTofDigitsCollection * m_besTofDigitsCollection
static NTuple::Item< double > m_nHits
static NTuple::Item< double > m_partIdMPV
static NTuple::Item< double > m_time1st1
static NTuple::Item< double > m_nDigiOut
BesTofHitsCollection * m_THC
static NTuple::Item< double > m_totalPhot1
static NTuple::Tuple * m_tupleTof1
static NTuple::Item< double > m_edepMPV
static NTuple::Item< double > m_time1st0
double GetBeamTime()
Definition: G4Svc.h:93
bool TofRootFlag()
Definition: G4Svc.h:126
G4int GetPartId()
Definition: ScintSingle.hh:44
vector< G4int > * GetHitIndexes()
Definition: ScintSingle.hh:47
G4int GetScinNb()
Definition: ScintSingle.hh:45
G4double GetEdep()
Definition: ScintSingle.hh:46
#define ns(x)
Definition: xmltok.c:1504

Referenced by BesTofDigitizer::Digitize().

◆ DirectPh()

void BesTofDigitizerBrV2::DirectPh ( G4int  partId,
G4ThreeVector  emtPos,
G4double &  pathL,
G4int &  forb 
)

Definition at line 463 of file BesTofDigitizerBrV2.cc.

464{
465 //std::cout << "BesTofDigitizerBrV2::DirectPh()" << std::endl;
466 const static G4double silicon_oil_index = 1.43;
467 const static G4double glass_index = 1.532;
468 //generation photon have random direction
469 //optical parameters of scintillation simulation
470 G4double cos_span=1-cos( asin(silicon_oil_index/m_refIndex) );
471 //G4double cos_spanEc = 1;
472 G4double dcos, ran;
473 ran=G4UniformRand();
474 //assuming uniform phi distribution, simulate cos distr only
475 dcos=cos_span*(ran*2.0-1.0);
476 //G4double dcosEc = cos_spanEc*(ran*2.0-1.0);
477 G4double r2=sqrt(emtPos.x()*emtPos.x()+emtPos.y()*emtPos.y());
478 G4int nRef=0;
479 G4double costheta,sintheta;
480 G4double theta,thetaR; // thetaR is scin to air full ref angle. about 39.265 degree.
481 costheta=dcos>0?(1-dcos):(1+dcos);
482 theta=acos(costheta);
483 sintheta=sin(theta);
484 thetaR=asin(1/m_refIndex);
485 G4double R1;
486 R1=Reflectivity(m_refIndex,1.0,theta);
487 G4double ratio1Mean=(CLHEP::pi)*25.5*25.5/(57.1+60.7)/25.0; //0.693657
488 G4double ratio2Mean=(19.5/25.5)*(19.5/25.5); //0.584775
489 G4double ratio1=G4RandGauss::shoot(ratio1Mean,0.015);
490 G4double ratio2=G4RandGauss::shoot(ratio2Mean,0.015);
491 G4double ratio3Mean=0.945*ratio1Mean*ratio2Mean;
492 G4double ratio3=G4RandGauss::shoot(ratio3Mean,0.015);
493
494 G4double R2=1-Reflectivity(m_refIndex,silicon_oil_index, glass_index, theta);
495 if (dcos > 0 && dcos != 1)
496 {
497 if (theta < thetaR) // theta < 39.265 degree
498 {
499 if (G4UniformRand()<ratio1) //coup1
500 {
501 if (G4UniformRand()<R2)
502 {
503 if (G4UniformRand()<ratio2) //PMT 39mm
504 {
505 forb=0;
506 pathL=(m_scinLength/2-emtPos.z())/costheta;
507 }
508 }
509 else
510 {
511 if (G4UniformRand()<ratio3)
512 {
513 forb=1;
514 pathL=(1.5*m_scinLength-emtPos.z())/costheta;
515 }
516 }
517 }
518 else //Air
519 {
520 if (G4UniformRand()<R1)
521 {
522 G4double tempran=G4UniformRand();
523 if (tempran<ratio3)
524 {
525 forb=1;
526 pathL=(1.5*m_scinLength-emtPos.z())/costheta;
527 }
528 else if (tempran>ratio1&&G4UniformRand()<R1&&G4UniformRand()<ratio3)
529 {
530 forb=0;
531 pathL=(2.5*m_scinLength-emtPos.z())/costheta;
532 }
533 }
534 }
535 }
536 else // 39.265 <= theta < 64.832
537 {
538 if (G4UniformRand()<ratio1) //coup1
539 {
540 if (G4UniformRand()<R2)
541 {
542 if (G4UniformRand()<ratio2) //PMT 39mm
543 {
544 forb=0;
545 pathL=(m_scinLength/2-emtPos.z())/costheta;
546 }
547 }
548 else
549 {
550 if (G4UniformRand()<ratio3)
551 {
552 forb=1;
553 pathL=(1.5*m_scinLength-emtPos.z())/costheta;
554 }
555 }
556 }
557 else //Air
558 {
559 G4double tempran=G4UniformRand();
560 if (tempran<ratio3)
561 {
562 forb=1;
563 pathL=(1.5*m_scinLength-emtPos.z())/costheta;
564 }
565 else if (tempran>ratio1 && G4UniformRand()<ratio3)
566 {
567 forb=0;
568 pathL=(2.5*m_scinLength-emtPos.z())/costheta;
569 }
570 }
571 }
572 }
573 else if (dcos < 0 && dcos != -1)
574 {
575 if (theta < thetaR) // theta < 39.265 degree
576 {
577 if (G4UniformRand()<ratio1) //coup1
578 {
579 if (G4UniformRand()<R2)
580 {
581 if (G4UniformRand()<ratio2) //PMT 39mm
582 {
583 forb=1;
584 pathL=(m_scinLength/2+emtPos.z())/costheta;
585 }
586 }
587 else
588 {
589 if (G4UniformRand()<ratio3)
590 {
591 forb=0;
592 pathL=(1.5*m_scinLength+emtPos.z())/costheta;
593 }
594 }
595 }
596 else //Air
597 {
598 if (G4UniformRand()<R1)
599 {
600 G4double tempran=G4UniformRand();
601 if (tempran<ratio3)
602 {
603 forb=0;
604 pathL=(1.5*m_scinLength+emtPos.z())/costheta;
605 }
606 else if (tempran>ratio1&&G4UniformRand()<R1&&G4UniformRand()<ratio3)
607 {
608 forb=1;
609 pathL=(2.5*m_scinLength+emtPos.z())/costheta;
610 }
611 }
612 }
613 }
614 else // 39.265 <= theta < 64.832
615 {
616 if (G4UniformRand()<ratio1) //coup1
617 {
618 if (G4UniformRand()<R2)
619 {
620 if (G4UniformRand()<ratio2) //PMT 39mm
621 {
622 forb=1;
623 pathL=(m_scinLength/2+emtPos.z())/costheta;
624 }
625 }
626 else
627 {
628 if (G4UniformRand()<ratio3)
629 {
630 forb=0;
631 pathL=(1.5*m_scinLength+emtPos.z())/costheta;
632 }
633 }
634 }
635 else //Air
636 {
637 G4double tempran=G4UniformRand();
638 if (tempran<ratio3)
639 {
640 forb=0;
641 pathL=(1.5*m_scinLength+emtPos.z())/costheta;
642 }
643 else if (tempran>ratio1 && G4UniformRand()<ratio3)
644 {
645 forb=1;
646 pathL=(2.5*m_scinLength+emtPos.z())/costheta;
647 }
648 }
649 }
650 }
651
652 G4double convFactor=180./3.1415926;
653 if (theta>asin(1)-thetaR)
654 {
655 G4double alpha = G4UniformRand()*asin(1.);
656 G4int nRef1 = pathL*sintheta*cos(alpha)/50.0+0.5;
657 G4int nRef2 = pathL*sintheta*sin(alpha)/58.9+0.5;
658 G4double beta1=acos(sintheta*cos(alpha));
659 G4double beta2=acos(sintheta*sin(alpha));
660 beta2 -= 3.75/convFactor;
661 G4double R21,R22;
662 R21=Reflectivity(m_refIndex,1.0,beta1);
663 R22=Reflectivity(m_refIndex,1.0,beta2);
664 for (G4int i=0;i<nRef1;i++)
665 {
666 if (G4UniformRand()<(1-R21) && G4UniformRand()<0.15)
667 pathL=-9;
668 }
669 for (G4int i=0;i<nRef2;i++)
670 {
671 if (G4UniformRand()<(1-R22) && G4UniformRand()<0.15)
672 pathL=-9;
673 }
674 }
675}
double sin(const BesAngle a)
Definition: BesAngle.h:210
double cos(const BesAngle a)
Definition: BesAngle.h:213
const double alpha
G4double Reflectivity(G4double n1, G4double n2, G4double theta)
float costheta

Referenced by TofPmtAccum().

◆ ReadData()

void BesTofDigitizerBrV2::ReadData ( )

Definition at line 62 of file BesTofDigitizerBrV2.cc.

63{
65
66 m_scinLength = tofPara->GetBr1L();
67 m_tau1 = tofPara->GetTau1();
68 m_tau2 = tofPara->GetTau2();
69 m_tau3 = tofPara->GetTau3();
70 m_tauRatio = tofPara->GetTauRatio();
71 m_refIndex = tofPara->GetRefIndex();
72 m_phNConst = tofPara->GetPhNConst();
73 m_Cpe2pmt = tofPara->GetCpe2pmt();
74 m_rAngle = tofPara->GetRAngle();
75 m_QE = tofPara->GetQE();
76 m_CE = tofPara->GetCE();
77 m_peCorFac = tofPara->GetPeCorFac();
78
79 m_ttsMean = tofPara->GetTTSmean();
80 m_ttsSigma = tofPara->GetTTSsigma();
81 m_Ce = tofPara->GetCe();
82 m_LLthresh = tofPara->GetLLthresh();
83 m_HLthresh = tofPara->GetHLthresh();
84 m_preGain = tofPara->GetPreGain();
85 m_noiseSigma = tofPara->GetNoiseSigma();
86
87
88}
static BesTofGeoParameter * GetInstance()

Referenced by BesTofDigitizerBrV2().

◆ Reflectivity() [1/2]

G4double BesTofDigitizerBrV2::Reflectivity ( G4double  n1,
G4double  n2,
G4double  n3,
G4double  theta 
)

Definition at line 416 of file BesTofDigitizerBrV2.cc.

417{
418 G4double I1,I2,I3,rp1,rs1,rp2,rs2,Rp1,Rs1,Rp2,Rs2,Rp,Rs;
419 G4double R=1.0;
420 //n1=m_refIndex;
421 //n2=1.0;
422 I1=theta;
423 if(I1<asin(n2/n1))
424 {
425 I2=asin(sin(I1)*(n1/n2));
426 rp1=(n1/cos(I1)-n2/cos(I2))/(n1/cos(I1)+n2/cos(I2));
427 rs1=(n1*cos(I1)-n2*cos(I2))/(n1*cos(I1)+n2*cos(I2));
428 Rp1=rp1*rp1;
429 Rs1=rs1*rs1;
430
431 I3=asin(sin(I1)*(n1/n3));
432 rp2=(n2/cos(I2)-n3/cos(I3))/(n2/cos(I2)+n3/cos(I3));
433 rs2=(n2*cos(I2)-n3*cos(I3))/(n2*cos(I2)+n3*cos(I3));
434 Rp2=rp2*rp2;
435 Rs2=rs2*rs2;
436 Rp=(Rp1+Rp2-2*Rp1*Rp2)/(1-Rp1*Rp2);
437 Rs=(Rs1+Rs2-2*Rs1*Rs2)/(1-Rs1*Rs2);
438 R=(Rp+Rs)/2.;
439 }
440 return R;
441}
int n2
Definition: SD0Tag.cxx:55
int n1
Definition: SD0Tag.cxx:54
complex_t R(double Q2, double M2, double G, double Mp2, double Mm2)
Definition: TUtil.h:27

◆ Reflectivity() [2/2]

G4double BesTofDigitizerBrV2::Reflectivity ( G4double  n1,
G4double  n2,
G4double  theta 
)

Definition at line 444 of file BesTofDigitizerBrV2.cc.

445{
446 G4double I1,I2,rp,rs,Rp,Rs;
447 G4double R=1.0;
448 //n1=m_refIndex;
449 //n2=1.0;
450 I1=theta;
451 if (I1<asin(n2/n1))
452 {
453 I2=asin(sin(I1)*(n1/n2));
454 rp=(n1/cos(I1)-n2/cos(I2))/(n1/cos(I1)+n2/cos(I2));
455 rs=(n1*cos(I1)-n2*cos(I2))/(n1*cos(I1)+n2*cos(I2));
456 Rp=rp*rp;
457 Rs=rs*rs;
458 R=(Rp+Rs)/2.;
459 }
460 return R;
461}

Referenced by DirectPh().

◆ Scintillation()

G4double BesTofDigitizerBrV2::Scintillation ( G4int  partId)

Definition at line 699 of file BesTofDigitizerBrV2.cc.

700{
701 G4double tmp_tauRatio,tmp_tau1,tmp_tau2,tmp_tau3;
702 tmp_tauRatio = m_tauRatio;
703 tmp_tau1 = m_tau1;
704 tmp_tau2 = m_tau2;
705 tmp_tau3 = m_tau3;
706 G4double UniformR = tmp_tauRatio/(1+tmp_tauRatio);
707 G4double EmissionTime;
708 if (G4UniformRand()>UniformR) {
709 while (1) {
710 EmissionTime = -tmp_tau2*log( G4UniformRand() );
711 if (G4UniformRand()-exp(EmissionTime/tmp_tau2-EmissionTime/tmp_tau1)>1.E-8)
712 break;
713 }
714 }
715 else EmissionTime = -tmp_tau3*log( G4UniformRand() );
716 return EmissionTime;
717}
EvtComplex exp(const EvtComplex &c)
Definition: EvtComplex.hh:252

Referenced by TofPmtAccum().

◆ TofPmtAccum()

void BesTofDigitizerBrV2::TofPmtAccum ( BesTofHit hit,
G4int  scinNb 
)

Definition at line 248 of file BesTofDigitizerBrV2.cc.

249{
251 //std::cout << "BesTofDigitizerBrV2::TofPmtAccum()" << std::endl;
252 G4double cvelScint=c_light/m_refIndex/1.07;
253 //Get information of this step
254 G4ThreeVector pos=hit->GetPos();
255 G4int trackIndex = hit->GetTrackIndex();
256 G4int partId =hit->GetPartId();
257 G4double edep=hit->GetEdep();
258 G4double stepL=hit->GetStepL();
259 G4double deltaT=hit->GetDeltaT();
260 G4double timeFlight=hit->GetTime()-m_beamTime;
261 //std::cout << "timeFlight: " << timeFlight << std::endl;
262 if (timeFlight < m_globalTime)
263 {
264 m_globalTime = timeFlight;
265 m_trackIndex = trackIndex;
266 }
267
268 G4ThreeVector pDirection=hit->GetPDirection();
269 G4double nx=pDirection.x();
270 G4double ny=pDirection.y();
271 G4double nz=pDirection.z();
272
273 //phNConst=(Nph/MeV)*(QE)*(CE)*(1-cos(crit));
274 // =10000 * 0.2 * 0.6 * (1-cos(39.265))=270.931
275 //asin(1/1.58) = 39.265248
276 //only the light has theta=0---39.265 can go out to PMT, the probability is computed with solid angle
277 //solid angle = phi*(1-cos(theta)), phi=2pi
278
279 //Cpe2pmt=cathode area/scint area
280 //peCorFac=correction factor for eff. of available Npe
281
282 //G4double thetaProb = 1-sqrt(m_refIndex*m_refIndex-1)/m_refIndex;
283 G4double thetaProb=1-cos( asin(1.43/m_refIndex));
284 //G4double thetaProbEc = 1-1/m_refIndex;
285
286 //number of photons generated in this step
287 G4double nMean, nPhoton;
288 //std::cout << "0 BirksLaw(): " << std::endl;
289 nMean = m_phNConst*BirksLaw(hit);
290 G4int runId = m_RealizationSvc->getRunId();
291 if( ( runId>=-11396 && runId<=-8093 ) || ( runId>-80000 && runId<=-23463 ) ) {
292 nMean = 9000.0*BirksLaw(hit);
293 }
294 //std::cout << "1 BirksLaw(): " << std::endl;
295
296 if(nMean>10)
297 {
298 G4double resolutionScale=1.;
299 G4double sigma=resolutionScale*sqrt(nMean);
300 nPhoton=G4int(G4RandGauss::shoot(nMean,sigma)+0.5);
301 }
302 else
303 nPhoton=G4int(G4Poisson(nMean));
304 //G4cout<<"nPhoton:"<<nPhoton<<G4endl;
305
306 G4int NphStep;
307 if (partId==1)
308 NphStep=G4int(nPhoton*thetaProb*m_rAngle*m_QE*m_CE*m_peCorFac);
309 //else
310 //NphStep=G4int(edep*m_phNConstEc*m_Cpe2pmtEc*0.66*m_QEEc*m_CE*m_peCorFac);
311 //introduce poission distribution of Npe
312 //G4double Navr = NphStep;
313 //G4double NpePoiss = G4double(G4Poisson(Navr));
314 //G4double adcW_corr =5.0;
315 //NphStep=G4int( (NpePoiss - Navr)*adcW_corr + Navr );
316
317 if (m_G4Svc->TofRootFlag())
318 m_NphAllSteps += NphStep;
319 //std::cout << "m_G4Svc->TofRootFlag(): " << m_G4Svc->TofRootFlag() << std::endl;
320
321 G4double ddS, ddT;
322 if (NphStep>0)
323 {
324 for (G4int i=0;i<NphStep;i++)
325 {
326 //uniform scintilation in each step
327 ddS=stepL*G4UniformRand();
328 ddT=deltaT*G4UniformRand();
329 G4ThreeVector emtPos;
330 emtPos.setX(pos.x() + nx*ddS);
331 emtPos.setY(pos.y() + ny*ddS);
332 emtPos.setZ(pos.z() + nz*ddS);
333
334 //check scinillation light whether to hit the pmt or not
335 //forb=0/1 for forward(z>0, east) and backward(z<0, west)
336 G4double pathL=0;
337 G4int forb;
338 DirectPh(partId, emtPos, pathL, forb);
339
340
341 //check if photon can reach PMT or not, after attenuation
342
343 G4double ran = G4UniformRand();
344 //G4double attenL = tofPara->GetAtten(scinNb);
345 //G4double attenL = 10.*(m_tofCaliSvc->BAtten(scinNb))/0.75; // cm into mm
346 G4double attenL = m_tofSimSvc->BarAttenLength(scinNb);
347 attenL = 10.*attenL/0.75; // cm into mm
348
349 if (pathL>0 && exp(-pathL/attenL) > ran)
350 {
351 //propagation time in scintillator
352 G4double scinSwim=pathL/cvelScint;
353 //scintillation timing
354 //G4double scinTime=GenPhoton(partId);
355 G4double scinTime=Scintillation(partId);
356
357 //PMT transit time
358 G4double transitTime=TransitTime();
359 //sum up all time components
360 G4double endTime= timeFlight + ddT + scinSwim + scinTime + transitTime;
361 //std::cout << "endtime: " << endTime << std::endl;
362
363 if (m_G4Svc->TofRootFlag())
364 {
365 //m_forb = forb;
366 //m_timeFlight = timeFlight;
367 //m_ddT = ddT;
368 m_scinSwim = scinSwim;
369 //m_scinTime = scinTime;
370 //m_transitTime = transitTime;
371 m_endTime = endTime;
372 m_tupleTof3->write();
373 }
374 //store timing into binned buffer
375 AccuSignal(endTime, forb);
376
377 //update 1st and last timings here
378 if (m_t1st[forb]>endTime) m_t1st[forb] = endTime;
379 if (m_tLast[forb]<endTime) m_tLast[forb]= endTime;
380 }
381 }
382 }
383}
TTree * sigma
void AccuSignal(G4double, G4int)
G4double Scintillation(G4int)
void DirectPh(G4int, G4ThreeVector, G4double &, G4int &)
G4double BirksLaw(BesTofHit *hit)
static NTuple::Tuple * m_tupleTof3
static NTuple::Item< double > m_NphAllSteps
static NTuple::Item< double > m_scinSwim
static NTuple::Item< double > m_endTime
ITofSimSvc * m_tofSimSvc
G4double GetDeltaT()
Definition: BesTofHit.hh:75
G4ThreeVector GetPDirection()
Definition: BesTofHit.hh:76
G4double GetTime()
Definition: BesTofHit.hh:74
G4ThreeVector GetPos()
Definition: BesTofHit.hh:73
G4int GetPartId()
Definition: BesTofHit.hh:68
G4int GetTrackIndex()
Definition: BesTofHit.hh:66
virtual const double BarAttenLength(unsigned int id)=0

Referenced by Digitize().

◆ TofPmtInit()

void BesTofDigitizerBrV2::TofPmtInit ( )

Definition at line 212 of file BesTofDigitizerBrV2.cc.

213{
214 Initialize();
215
216 m_t1st[0]=100;
217 m_t1st[1]=100;
218 m_tLast[0]=0.;
219 m_tLast[1]=0;
220 m_totalPhot[0]=0;
221 m_totalPhot[1]=0;
222 for (G4int i=0;i<2;i++)
223 for (G4int j=0;j<m_profBinN;j++)
224 m_nPhot[j][i]=0;
225
226 if (m_G4Svc->TofRootFlag())
227 {
228 m_partId = -9;
229 m_scinNb = -9;
230 m_edep = 0;
231 m_nHits = 0;
232 m_time1st0 = 100;
233 m_time1st1 = 100;
234 m_timelast0 = 0;
235 m_timelast1 = 0;
236 m_totalPhot0 = 0;
237 m_totalPhot1 = 0;
238 m_NphAllSteps = 0;
239 m_max0 = 0;
240 m_max1 = 0;
241 m_tdc0 = -999;
242 m_adc0 = -999;
243 m_tdc1 = -999;
244 m_adc1 = -999;
245 }
246}
static NTuple::Item< double > m_tdc1
static NTuple::Item< double > m_max0
static NTuple::Item< double > m_tdc0
static NTuple::Item< double > m_adc1
static NTuple::Item< double > m_adc0
static NTuple::Item< double > m_max1

Referenced by Digitize().

◆ TofPmtRspns()

void BesTofDigitizerBrV2::TofPmtRspns ( G4int  scinNb)

Definition at line 763 of file BesTofDigitizerBrV2.cc.

764{
765 //std::cout << "BesTofDigitizerBrV2::TofPmtRspns()" << std::endl;
767 //to generate PMT signal shape for single photoelectron.
768 //only one time for a job.
769 G4double snpe[m_snpeBinN][2];
770
771 //Model: f(t)=Gain*mv_1pe* t**2 * exp-(t/tau)**2/normal-const
772 //normalization const =sqrt(pi)*tau*tau*tau/4.0
773 //G4double tau = m_riseTime;
774 G4double norma_const;
775 G4double echarge=1.6e-7; //in unit of PC
776
777 //time profile of pmt signals for Back and Forward PMT.
778 G4double profPmt[m_profBinN][2];
779
780 G4double t;
781 G4double tau;
782 G4int n1, n2, ii;
783 G4int phtn;
784 // forb = 0, east
785 for (G4int i=0;i<m_snpeBinN;i++)
786 {
787 tau = tofPara->GetBrERiseTime(scinNb);
788 norma_const = sqrt(CLHEP::pi)*tau*tau*tau/4.0; //in unit of ns**3
789 t = (i+1)*m_timeBinSize;
790 snpe[i][0] = m_Ce*t*t*exp(- (t/tau) * (t/tau) )/norma_const;// Pulse of one single photoelectron
791 }
792 // forb = 1, west
793 for (G4int i=0;i<m_snpeBinN;i++)
794 {
795 tau = tofPara->GetBrWRiseTime(scinNb);
796 norma_const = sqrt(CLHEP::pi)*tau*tau*tau/4.0; //in unit of ns**3
797 t = (i+1)*m_timeBinSize;
798 snpe[i][1] = m_Ce*t*t*exp(- (t/tau) * (t/tau) )/norma_const;
799 }
800 //for barrel and endcap tof: fb=2 or 1
801 G4int fb=2;
802 G4double Npoisson;
803 G4double pmtGain0,pmtGain,relativeGain,smearPMTgain;
804 G4double smearADC[2] = {0.,0.};
805 G4int runId = m_RealizationSvc->getRunId();
806 pmtGain0 = m_tofSimSvc->BarPMTGain();
807
808// if(runId>=-80000 && runId<=-9484)
809// {
810// // After TOF HV adjustment, PMT Gain of barrel TOF in MC is set to 5.E5
811// // High/Low Threshold for barrel: 500/125 mV
812// pmtGain0 = 6.E5;
813// }
814// else
815// {
816// // Before TOF HV adjustment, PMT Gain of barrel TOF in MC is set to 2.5E5
817// // High/Low Threshold for barrel: 600/150 mV
818// pmtGain0 = 5.E5;
819// }
820
821 G4double timeSmear = G4RandGauss::shoot(0,0.020);
822 if (runId>=-22913 && runId<=-20448) {//for 2011-psipp(20448-22913), smear barrel TOF resolution to ~78ps
823 timeSmear = G4RandGauss::shoot(0,0.040);
824 }
825 else if( ( runId>=-11396 && runId<=-8093 ) || ( runId>-80000 && runId<=-23463 ) ) {
826 timeSmear = G4RandGauss::shoot(0,0.025);
827 }
828
829 for (G4int j=0; j<fb; j++)
830 {
831 if (m_totalPhot[j] > 0)
832 {
833 n1=G4int(m_t1st[j]/m_timeBinSize);
834 n2=G4int(m_tLast[j]/m_timeBinSize);
835 //std::cout << "n1: " << n1 << std::endl;
836
837 for (G4int i=0;i<m_profBinN;i++)
838 profPmt[i][j]=0.0;
839
840 //generate PMT pulse
842 for (G4int i=n1;i<n2;i++)
843 {
844 phtn=m_nPhot[i][j];
845 if (phtn>0)
846 {
847 while(1) {
848 Npoisson = G4Poisson(10.0);
849 if(Npoisson>0.) break;
850 }
851 while(1) {
852 //pmtGain = j ? tofPara->GetBrWPMTgain(scinNb) : tofPara->GetBrEPMTgain(scinNb);
853 //relativeGain = j ? m_tofCaliSvc->BGainBackward(scinNb) : m_tofCaliSvc->BGainForward(scinNb);
854 relativeGain = j ? m_tofSimSvc->BarGain2(scinNb) : m_tofSimSvc->BarGain1(scinNb);
855 pmtGain = pmtGain0*relativeGain;
856 smearPMTgain = G4RandGauss::shoot(pmtGain,pmtGain/sqrt(Npoisson));
857 //smearPMTgain = pmtGain;
858 if(smearPMTgain>0) break;
859 }
860 smearADC[j] += phtn*smearPMTgain;
861
862 for (G4int ihst=0; ihst<m_snpeBinN; ihst++)
863 {
864 ii=i+ihst;
865 if (ii<m_profBinN)
866 profPmt[ii][j] += smearPMTgain*phtn*snpe[ihst][j];
867 else
868 break;
869 }
870 }
871 }
872
873 //add preamplifier and noise
874 for (int i=0;i<m_profBinN;i++)
875 {
876 if (profPmt[i][j]>0)
877 profPmt[i][j] = m_preGain*profPmt[i][j]+G4RandGauss::shoot(0,m_noiseSigma);
878 }
879
880 //get pulse height
881 G4double max=0;
882 for (int i=n1;i<m_profBinN;i++)
883 {
884 if (profPmt[i][j]>max)
885 max=profPmt[i][j];
886 }
887 if (m_G4Svc->TofRootFlag())
888 {
889 if (j==0) m_max0=max;
890 else m_max1=max;
891 }
892
893 G4double tmp_HLthresh, tmp_LLthresh, adcFactor;
894 G4double ratio;
895
896 if(runId>=-80000 && runId<=-9484)
897 {
898 // After TOF HV adjustment, PMT Gain of barrel TOF in MC is set to 5.E5
899 // High/Low Threshold for barrel: 500/125 mV
900 tmp_HLthresh = m_HLthresh;
901 tmp_LLthresh = m_LLthresh;
902 adcFactor = 5.89;
903 }
904 else
905 {
906 // Before TOF HV adjustment, PMT Gain of barrel TOF in MC is set to 2.5E5
907 // High/Low Threshold for barrel: 600/150 mV
908 tmp_HLthresh = 600;
909 tmp_LLthresh = 150;
910 adcFactor = 4.8;
911 }
912// if(runId>=-80000 && runId<=-9484)
913// {
914// ratio=2.16*1923.8/2197.8;
915// }
916// else
917// {
918// ratio = 2.11*2437.0/3102.9;
919// }
920 tmp_HLthresh = m_tofSimSvc->BarHighThres();
921 tmp_LLthresh = m_tofSimSvc->BarLowThres();
922 ratio = m_tofSimSvc->BarConstant();
923
924 //get final tdc and adc
925 if (max>=tmp_HLthresh)
926 {
927 for (int i=0;i<m_profBinN;i++)
928 {
929 if ( profPmt[i][j] >= tmp_LLthresh )
930 {
931 m_TDC[j] = i*m_timeBinSize + timeSmear +G4RandGauss::shoot(0,0.025); // Adding Electronical Uncertainty of 25ps
932// hhliu 20140724, setting tdc-smear for inner and outer separately
933 if( ( runId>=-11396 && runId<=-8093 ) || ( runId>-80000 && runId<=-23463 ) ) {
934 if( scinNb<88 ) {
935 m_TDC[j] = i*m_timeBinSize + timeSmear +G4RandGauss::shoot(0,0.055);
936 }
937 else {
938 m_TDC[j] = i*m_timeBinSize + timeSmear +G4RandGauss::shoot(0,0.062);
939 }
940 }
941
942 if( m_G4Svc->TofSaturationFlag())
943 {
944 //get ADC[j] using tofQElecSvc
945 double x = m_preGain*smearADC[j]*echarge*ratio;
946 if (j==0)
947 {
948 m_ADC[j] = m_tofQElecSvc->BQChannel1(scinNb,x);
949 }
950 else
951 {
952 m_ADC[j] = m_tofQElecSvc->BQChannel2(scinNb,x);
953 }
954 }
955 else
956 m_ADC[j] = m_preGain*smearADC[j]*echarge*adcFactor;
957
958
959 if (m_G4Svc->TofRootFlag())
960 {
961 if (j==0) {
962 m_tdc0 = m_TDC[0];
963 m_adc0 = m_ADC[0];
964 }
965 else {
966 m_tdc1 = m_TDC[1];
967 m_adc1 = m_ADC[1];
968 }
969 }
970 break;
971 }
972 }
973 }
974 }
975 }
976}
const G4int m_snpeBinN
Double_t x[10]
TTree * t
Definition: binning.cxx:23
ITofQElecSvc * m_tofQElecSvc
G4double GetBrWRiseTime(int scinNb)
G4double GetBrERiseTime(int scinNb)
bool TofSaturationFlag()
Definition: G4Svc.h:130
virtual const double BQChannel2(int id, double q)=0
virtual const double BQChannel1(int id, double q)=0
virtual const double BarConstant()=0
virtual const double BarLowThres()=0
virtual const double BarPMTGain()=0
virtual const double BarGain2(unsigned int id)=0
virtual const double BarHighThres()=0
virtual const double BarGain1(unsigned int id)=0

Referenced by Digitize().

◆ TransitTime()

G4double BesTofDigitizerBrV2::TransitTime ( )

Definition at line 746 of file BesTofDigitizerBrV2.cc.

747{
748 //get time transit spread
749 return G4RandGauss::shoot(m_ttsMean,m_ttsSigma);
750}

Referenced by TofPmtAccum().


The documentation for this class was generated from the following files: