BOSS 7.0.3
BESIII Offline Software System
Loading...
Searching...
No Matches
BesTofDigitizerBrV2.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: BesTofDigitizerBrV2.cc
11
12#include "BesTofDigitizerBrV2.hh"
13#include "BesTofHit.hh"
14#include "G4DigiManager.hh"
15#include "G4RunManager.hh"
16#include "Randomize.hh"
17#include "G4Poisson.hh"
18#include "BesTofDigi.hh"
19#include "BesTofGeoParameter.hh"
20#include "GaudiKernel/ISvcLocator.h"
21#include "GaudiKernel/Bootstrap.h"
22#include "GaudiKernel/IDataProviderSvc.h"
23#include "G4Svc/IG4Svc.h"
24#include "G4Svc/G4Svc.h"
25#include "TMath.h"
26
28{
29 ReadData();
30 m_timeBinSize=0.005;
31
32 //retrieve G4Svc
33 ISvcLocator* svcLocator = Gaudi::svcLocator();
34 IG4Svc* tmpSvc;
35 StatusCode sc = svcLocator->service("G4Svc", tmpSvc);
36 if(!sc.isSuccess())
37 {
38 std::cout << " Could not initialize Realization Service in BesTofDigitizerBrV2" << std::endl;
39 }
40 else
41 {
42 m_G4Svc = dynamic_cast<G4Svc *>(tmpSvc);
43 }
44
45 //retrieve RealizationSvc
46 IRealizationSvc *tmpReal;
47 StatusCode scReal = svcLocator->service("RealizationSvc",tmpReal);
48 if (!scReal.isSuccess())
49 {
50 std::cout << " Could not initialize Realization Service in BesTofDigitizerBrV2" << std::endl;
51 }
52 else
53 {
54 m_RealizationSvc = dynamic_cast<RealizationSvc*>(tmpReal);
55 }
56
57
58}
59
61{
63
64 m_scinLength = tofPara->GetBr1L();
65 m_tau1 = tofPara->GetTau1();
66 m_tau2 = tofPara->GetTau2();
67 m_tau3 = tofPara->GetTau3();
68 m_tauRatio = tofPara->GetTauRatio();
69 m_refIndex = tofPara->GetRefIndex();
70 m_phNConst = tofPara->GetPhNConst();
71 m_Cpe2pmt = tofPara->GetCpe2pmt();
72 m_rAngle = tofPara->GetRAngle();
73 m_QE = tofPara->GetQE();
74 m_CE = tofPara->GetCE();
75 m_peCorFac = tofPara->GetPeCorFac();
76
77 m_ttsMean = tofPara->GetTTSmean();
78 m_ttsSigma = tofPara->GetTTSsigma();
79 m_Ce = tofPara->GetCe();
80 m_LLthresh = tofPara->GetLLthresh();
81 m_HLthresh = tofPara->GetHLthresh();
82 m_preGain = tofPara->GetPreGain();
83 m_noiseSigma = tofPara->GetNoiseSigma();
84
85
86}
87
89{
90 ;
91}
92
94{
95 m_beamTime = m_G4Svc->GetBeamTime() * ns;
97
98 G4DigiManager* digiManager = G4DigiManager::GetDMpointer();
99 G4int THCID = digiManager->GetHitsCollectionID("BesTofHitsCollection");
100 m_THC = (BesTofHitsCollection*) (digiManager->GetHitsCollection(THCID));
101
102 if (m_G4Svc->TofRootFlag())
103 {
104 m_eTotal = 0;
105 m_nDigi = 0;
106 m_partIdMPV = -9;
107 m_scinNbMPV = -9;
108 m_edepMPV = 0;
109 m_nDigiOut = 0;
110 }
111
112 if (m_THC)
113 {
114 //for each digi, compute TDC and ADC
115 G4int partId, scinNb, nHits;
116 G4double edep;
117 BesTofHit* hit;
118 partId=scint->GetPartId();
119 scinNb=scint->GetScinNb();
120 edep = scint->GetEdep();
121 nHits=scint->GetHitIndexes()->size();
122
123
124 //std::cout << "BesTofDigitizerBrV2 Partid scinNb " << partId << " " << scinNb << std::endl;
125 //cout << "*** scinNb:" << scinNb << " *** " << m_tofCaliSvc->BAtten(scinNb) << "***" << endl;
126 //cout << "*****scinNb:"<< scinNb << " ***** A2:" << m_tofCaliSvc->BGainBackward(scinNb)
127 // << " ***** A1:" << m_tofCaliSvc->BGainForward(scinNb) << endl;
128
129 TofPmtInit();
130
131 //fill tof Ntuple
132 if (m_G4Svc->TofRootFlag())
133 {
134 if (edep>m_edepMPV)
135 {
136 m_partIdMPV = partId;
137 m_scinNbMPV = scinNb;
138 m_edepMPV = edep;
139 }
140 m_eTotal += edep;
141 m_nDigi ++;
142
143 m_partId = partId;
144 m_scinNb = scinNb;
145 m_edep = edep;
146 m_nHits = nHits;
147 }
148
149 if (edep>0.01)
150 {
151 for (G4int j=0;j<nHits;j++)
152 {
153 hit= (*m_THC)[( *(scint->GetHitIndexes()) )[j]];
154 TofPmtAccum(hit, scinNb);
155 }
156 if (m_G4Svc->TofRootFlag())
157 {
158 m_time1st0=m_t1st[0];
159 m_time1st1=m_t1st[1];
160 m_timelast0=m_tLast[0];
161 m_timelast1=m_tLast[1];
162 m_totalPhot0=m_totalPhot[0];
163 m_totalPhot1=m_totalPhot[1];
164 }
165 //get final tdc and adc
166 TofPmtRspns(scinNb);
167 //G4cout<<"pre-cut " << partId << "\nadc0:"<<m_ADC[0]<<"; adc1:"
168 // <<m_ADC[1]<<"; tdc0:"<<m_TDC[0]<<"; tdc1:"
169 // <<m_TDC[1]<<G4endl;
170
171 G4double temp0 = m_ADC[0]+m_TDC[0];
172 G4double temp1 = m_ADC[1]+m_TDC[1];
173 //cout << "partid: " << partId << " temp0: " << temp0 << " temp1: " << temp1 << endl;
174 //if ( partId==1 && m_ADC[0]>255 && m_ADC[1]>255 && m_TDC[0]>0. && m_TDC[1]>0.)
175 if ( partId==1 && (temp0 > -1998. || temp1 > -1998.))
176 {
177 //const double MAX_ADC = 8191 * 0.3; // channel set up to 8192 will lead to overflow.
178 BesTofDigi* digi = new BesTofDigi;
180 digi->SetPartId(partId);
181 digi->SetScinNb(scinNb);
182 digi->SetForwADC( m_ADC[0]) ;
183 digi->SetBackADC( m_ADC[1]) ;
184 if (m_TDC[0]>0)
185 m_TDC[0] = m_TDC[0]+m_beamTime;
186 if (m_TDC[1]>0)
187 m_TDC[1] = m_TDC[1]+m_beamTime;
188 digi->SetForwTDC( m_TDC[0]) ;
189 digi->SetBackTDC( m_TDC[1]) ;
190 //G4cout<<"+++++++++++++++++++++++++++++barrel\nadc0:"<<m_ADC[0]<<"; adc1:"
191 // <<m_ADC[1]<<"; tdc0:"<<m_TDC[0]<<"; tdc1:"
192 // <<m_TDC[1]<<G4endl;
193
194 m_besTofDigitsCollection->insert(digi);
195
196 if (m_G4Svc->TofRootFlag())
197 m_nDigiOut++;
198
199 }
200 if (m_G4Svc->TofRootFlag())
201 m_tupleTof1->write();
202
203 }
204
205 if (m_G4Svc->TofRootFlag())
206 m_tupleTof2->write();
207 }
208}
209
211{
212 Initialize();
213
214 m_t1st[0]=100;
215 m_t1st[1]=100;
216 m_tLast[0]=0.;
217 m_tLast[1]=0;
218 m_totalPhot[0]=0;
219 m_totalPhot[1]=0;
220 for (G4int i=0;i<2;i++)
221 for (G4int j=0;j<m_profBinN;j++)
222 m_nPhot[j][i]=0;
223
224 if (m_G4Svc->TofRootFlag())
225 {
226 m_partId = -9;
227 m_scinNb = -9;
228 m_edep = 0;
229 m_nHits = 0;
230 m_time1st0 = 100;
231 m_time1st1 = 100;
232 m_timelast0 = 0;
233 m_timelast1 = 0;
234 m_totalPhot0 = 0;
235 m_totalPhot1 = 0;
236 m_NphAllSteps = 0;
237 m_max0 = 0;
238 m_max1 = 0;
239 m_tdc0 = -999;
240 m_adc0 = -999;
241 m_tdc1 = -999;
242 m_adc1 = -999;
243 }
244}
245
247{
249 //std::cout << "BesTofDigitizerBrV2::TofPmtAccum()" << std::endl;
250 G4double cvelScint=c_light/m_refIndex/1.07;
251 //Get information of this step
252 G4ThreeVector pos=hit->GetPos();
253 G4int trackIndex = hit->GetTrackIndex();
254 G4int partId =hit->GetPartId();
255 G4double edep=hit->GetEdep();
256 G4double stepL=hit->GetStepL();
257 G4double deltaT=hit->GetDeltaT();
258 G4double timeFlight=hit->GetTime()-m_beamTime;
259 //std::cout << "timeFlight: " << timeFlight << std::endl;
260 if (timeFlight < m_globalTime)
261 {
262 m_globalTime = timeFlight;
263 m_trackIndex = trackIndex;
264 }
265
266 G4ThreeVector pDirection=hit->GetPDirection();
267 G4double nx=pDirection.x();
268 G4double ny=pDirection.y();
269 G4double nz=pDirection.z();
270
271 //phNConst=(Nph/MeV)*(QE)*(CE)*(1-cos(crit));
272 // =10000 * 0.2 * 0.6 * (1-cos(39.265))=270.931
273 //asin(1/1.58) = 39.265248
274 //only the light has theta=0---39.265 can go out to PMT, the probability is computed with solid angle
275 //solid angle = phi*(1-cos(theta)), phi=2pi
276
277 //Cpe2pmt=cathode area/scint area
278 //peCorFac=correction factor for eff. of available Npe
279
280 //G4double thetaProb = 1-sqrt(m_refIndex*m_refIndex-1)/m_refIndex;
281 G4double thetaProb=1-cos( asin(1.43/m_refIndex));
282 //G4double thetaProbEc = 1-1/m_refIndex;
283
284 //number of photons generated in this step
285 G4double nMean, nPhoton;
286 //std::cout << "0 BirksLaw(): " << std::endl;
287 nMean = m_phNConst*BirksLaw(hit);
288 G4int runId = m_RealizationSvc->getRunId();
289 if( ( runId>=-11396 && runId<=-8093 ) || ( runId>-80000 && runId<=-23463 ) ) {
290 nMean = 9000.0*BirksLaw(hit);
291 }
292 //std::cout << "1 BirksLaw(): " << std::endl;
293
294 if(nMean>10)
295 {
296 G4double resolutionScale=1.;
297 G4double sigma=resolutionScale*sqrt(nMean);
298 nPhoton=G4int(G4RandGauss::shoot(nMean,sigma)+0.5);
299 }
300 else
301 nPhoton=G4int(G4Poisson(nMean));
302 //G4cout<<"nPhoton:"<<nPhoton<<G4endl;
303
304 G4int NphStep;
305 if (partId==1)
306 NphStep=G4int(nPhoton*thetaProb*m_rAngle*m_QE*m_CE*m_peCorFac);
307 //else
308 //NphStep=G4int(edep*m_phNConstEc*m_Cpe2pmtEc*0.66*m_QEEc*m_CE*m_peCorFac);
309 //introduce poission distribution of Npe
310 //G4double Navr = NphStep;
311 //G4double NpePoiss = G4double(G4Poisson(Navr));
312 //G4double adcW_corr =5.0;
313 //NphStep=G4int( (NpePoiss - Navr)*adcW_corr + Navr );
314
315 if (m_G4Svc->TofRootFlag())
316 m_NphAllSteps += NphStep;
317 //std::cout << "m_G4Svc->TofRootFlag(): " << m_G4Svc->TofRootFlag() << std::endl;
318
319 G4double ddS, ddT;
320 if (NphStep>0)
321 {
322 for (G4int i=0;i<NphStep;i++)
323 {
324 //uniform scintilation in each step
325 ddS=stepL*G4UniformRand();
326 ddT=deltaT*G4UniformRand();
327 G4ThreeVector emtPos;
328 emtPos.setX(pos.x() + nx*ddS);
329 emtPos.setY(pos.y() + ny*ddS);
330 emtPos.setZ(pos.z() + nz*ddS);
331
332 //check scinillation light whether to hit the pmt or not
333 //forb=0/1 for forward(z>0, east) and backward(z<0, west)
334 G4double pathL=0;
335 G4int forb;
336 DirectPh(partId, emtPos, pathL, forb);
337
338
339 //check if photon can reach PMT or not, after attenuation
340
341 G4double ran = G4UniformRand();
342 //G4double attenL = tofPara->GetAtten(scinNb);
343 //G4double attenL = 10.*(m_tofCaliSvc->BAtten(scinNb))/0.75; // cm into mm
344 G4double attenL = m_tofSimSvc->BarAttenLength(scinNb);
345 attenL = 10.*attenL/0.75; // cm into mm
346
347 if (pathL>0 && exp(-pathL/attenL) > ran)
348 {
349 //propagation time in scintillator
350 G4double scinSwim=pathL/cvelScint;
351 //scintillation timing
352 //G4double scinTime=GenPhoton(partId);
353 G4double scinTime=Scintillation(partId);
354
355 //PMT transit time
356 G4double transitTime=TransitTime();
357 //sum up all time components
358 G4double endTime= timeFlight + ddT + scinSwim + scinTime + transitTime;
359 //std::cout << "endtime: " << endTime << std::endl;
360
361 if (m_G4Svc->TofRootFlag())
362 {
363 //m_forb = forb;
364 //m_timeFlight = timeFlight;
365 //m_ddT = ddT;
366 m_scinSwim = scinSwim;
367 //m_scinTime = scinTime;
368 //m_transitTime = transitTime;
369 m_endTime = endTime;
370 m_tupleTof3->write();
371 }
372 //store timing into binned buffer
373 AccuSignal(endTime, forb);
374
375 //update 1st and last timings here
376 if (m_t1st[forb]>endTime) m_t1st[forb] = endTime;
377 if (m_tLast[forb]<endTime) m_tLast[forb]= endTime;
378 }
379 }
380 }
381}
382
384{
385 const G4double kappa = 0.015*cm/MeV;
386 const G4String brMaterial = "BC408";
387 G4double dE = hit->GetEdep();
388 //G4cout << "The edep is "<< dE << G4endl;
389 G4double dX = hit->GetStepL();
390 //G4Material* materiral = hit->GetMaterial();
391 G4double charge = hit->GetCharge();
392 G4double cor_dE = dE;
393 //if((materiral->GetName()==brMaterial) && charge!=0.&& dX!=0.)
394 if(charge!=0.&& dX!=0.)
395 {
396 cor_dE = dE/(1+kappa*dE/dX);
397 //if(dE>20)
398 //{
399 // G4cout << "\n dE > 20. Details are below:" << G4endl;
400 // G4cout << "dE/dx:" << dE/dX << G4endl;
401 // G4cout << "dE:" << dE << "; dX:" << dX << G4endl;
402 // G4cout << "It is BC408. cor_dE is " << cor_dE << G4endl;
403 // G4double ratio = cor_dE/dE;
404 // G4cout << "The ratio cor_dE/edep is "<< ratio << G4endl;
405 //}
406 //G4cout << "It is BC408. cor_dE is " << cor_dE << G4endl;
407 //G4double ratio = cor_dE/dE;
408 //G4cout << "The ratio cor_dE/edep is "<< ratio << G4endl;
409 }
410 return cor_dE;
411
412}
413
414G4double BesTofDigitizerBrV2::Reflectivity(G4double n1,G4double n2,G4double n3,G4double theta)
415{
416 G4double I1,I2,I3,rp1,rs1,rp2,rs2,Rp1,Rs1,Rp2,Rs2,Rp,Rs;
417 G4double R=1.0;
418 //n1=m_refIndex;
419 //n2=1.0;
420 I1=theta;
421 if(I1<asin(n2/n1))
422 {
423 I2=asin(sin(I1)*(n1/n2));
424 rp1=(n1/cos(I1)-n2/cos(I2))/(n1/cos(I1)+n2/cos(I2));
425 rs1=(n1*cos(I1)-n2*cos(I2))/(n1*cos(I1)+n2*cos(I2));
426 Rp1=rp1*rp1;
427 Rs1=rs1*rs1;
428
429 I3=asin(sin(I1)*(n1/n3));
430 rp2=(n2/cos(I2)-n3/cos(I3))/(n2/cos(I2)+n3/cos(I3));
431 rs2=(n2*cos(I2)-n3*cos(I3))/(n2*cos(I2)+n3*cos(I3));
432 Rp2=rp2*rp2;
433 Rs2=rs2*rs2;
434 Rp=(Rp1+Rp2-2*Rp1*Rp2)/(1-Rp1*Rp2);
435 Rs=(Rs1+Rs2-2*Rs1*Rs2)/(1-Rs1*Rs2);
436 R=(Rp+Rs)/2.;
437 }
438 return R;
439}
440
441
442G4double BesTofDigitizerBrV2::Reflectivity(G4double n1,G4double n2,G4double theta)
443{
444 G4double I1,I2,rp,rs,Rp,Rs;
445 G4double R=1.0;
446 //n1=m_refIndex;
447 //n2=1.0;
448 I1=theta;
449 if (I1<asin(n2/n1))
450 {
451 I2=asin(sin(I1)*(n1/n2));
452 rp=(n1/cos(I1)-n2/cos(I2))/(n1/cos(I1)+n2/cos(I2));
453 rs=(n1*cos(I1)-n2*cos(I2))/(n1*cos(I1)+n2*cos(I2));
454 Rp=rp*rp;
455 Rs=rs*rs;
456 R=(Rp+Rs)/2.;
457 }
458 return R;
459}
460
461void BesTofDigitizerBrV2::DirectPh(G4int partId, G4ThreeVector emtPos, G4double& pathL, G4int& forb)
462{
463 //std::cout << "BesTofDigitizerBrV2::DirectPh()" << std::endl;
464 const static G4double silicon_oil_index = 1.43;
465 const static G4double glass_index = 1.532;
466 //generation photon have random direction
467 //optical parameters of scintillation simulation
468 G4double cos_span=1-cos( asin(silicon_oil_index/m_refIndex) );
469 //G4double cos_spanEc = 1;
470 G4double dcos, ran;
471 ran=G4UniformRand();
472 //assuming uniform phi distribution, simulate cos distr only
473 dcos=cos_span*(ran*2.0-1.0);
474 //G4double dcosEc = cos_spanEc*(ran*2.0-1.0);
475 G4double r2=sqrt(emtPos.x()*emtPos.x()+emtPos.y()*emtPos.y());
476 G4int nRef=0;
477 G4double costheta,sintheta;
478 G4double theta,thetaR; // thetaR is scin to air full ref angle. about 39.265 degree.
479 costheta=dcos>0?(1-dcos):(1+dcos);
480 theta=acos(costheta);
481 sintheta=sin(theta);
482 thetaR=asin(1/m_refIndex);
483 G4double R1;
484 R1=Reflectivity(m_refIndex,1.0,theta);
485 G4double ratio1Mean=(CLHEP::pi)*25.5*25.5/(57.1+60.7)/25.0; //0.693657
486 G4double ratio2Mean=(19.5/25.5)*(19.5/25.5); //0.584775
487 G4double ratio1=G4RandGauss::shoot(ratio1Mean,0.015);
488 G4double ratio2=G4RandGauss::shoot(ratio2Mean,0.015);
489 G4double ratio3Mean=0.945*ratio1Mean*ratio2Mean;
490 G4double ratio3=G4RandGauss::shoot(ratio3Mean,0.015);
491
492 G4double R2=1-Reflectivity(m_refIndex,silicon_oil_index, glass_index, theta);
493 if (dcos > 0 && dcos != 1)
494 {
495 if (theta < thetaR) // theta < 39.265 degree
496 {
497 if (G4UniformRand()<ratio1) //coup1
498 {
499 if (G4UniformRand()<R2)
500 {
501 if (G4UniformRand()<ratio2) //PMT 39mm
502 {
503 forb=0;
504 pathL=(m_scinLength/2-emtPos.z())/costheta;
505 }
506 }
507 else
508 {
509 if (G4UniformRand()<ratio3)
510 {
511 forb=1;
512 pathL=(1.5*m_scinLength-emtPos.z())/costheta;
513 }
514 }
515 }
516 else //Air
517 {
518 if (G4UniformRand()<R1)
519 {
520 G4double tempran=G4UniformRand();
521 if (tempran<ratio3)
522 {
523 forb=1;
524 pathL=(1.5*m_scinLength-emtPos.z())/costheta;
525 }
526 else if (tempran>ratio1&&G4UniformRand()<R1&&G4UniformRand()<ratio3)
527 {
528 forb=0;
529 pathL=(2.5*m_scinLength-emtPos.z())/costheta;
530 }
531 }
532 }
533 }
534 else // 39.265 <= theta < 64.832
535 {
536 if (G4UniformRand()<ratio1) //coup1
537 {
538 if (G4UniformRand()<R2)
539 {
540 if (G4UniformRand()<ratio2) //PMT 39mm
541 {
542 forb=0;
543 pathL=(m_scinLength/2-emtPos.z())/costheta;
544 }
545 }
546 else
547 {
548 if (G4UniformRand()<ratio3)
549 {
550 forb=1;
551 pathL=(1.5*m_scinLength-emtPos.z())/costheta;
552 }
553 }
554 }
555 else //Air
556 {
557 G4double tempran=G4UniformRand();
558 if (tempran<ratio3)
559 {
560 forb=1;
561 pathL=(1.5*m_scinLength-emtPos.z())/costheta;
562 }
563 else if (tempran>ratio1 && G4UniformRand()<ratio3)
564 {
565 forb=0;
566 pathL=(2.5*m_scinLength-emtPos.z())/costheta;
567 }
568 }
569 }
570 }
571 else if (dcos < 0 && dcos != -1)
572 {
573 if (theta < thetaR) // theta < 39.265 degree
574 {
575 if (G4UniformRand()<ratio1) //coup1
576 {
577 if (G4UniformRand()<R2)
578 {
579 if (G4UniformRand()<ratio2) //PMT 39mm
580 {
581 forb=1;
582 pathL=(m_scinLength/2+emtPos.z())/costheta;
583 }
584 }
585 else
586 {
587 if (G4UniformRand()<ratio3)
588 {
589 forb=0;
590 pathL=(1.5*m_scinLength+emtPos.z())/costheta;
591 }
592 }
593 }
594 else //Air
595 {
596 if (G4UniformRand()<R1)
597 {
598 G4double tempran=G4UniformRand();
599 if (tempran<ratio3)
600 {
601 forb=0;
602 pathL=(1.5*m_scinLength+emtPos.z())/costheta;
603 }
604 else if (tempran>ratio1&&G4UniformRand()<R1&&G4UniformRand()<ratio3)
605 {
606 forb=1;
607 pathL=(2.5*m_scinLength+emtPos.z())/costheta;
608 }
609 }
610 }
611 }
612 else // 39.265 <= theta < 64.832
613 {
614 if (G4UniformRand()<ratio1) //coup1
615 {
616 if (G4UniformRand()<R2)
617 {
618 if (G4UniformRand()<ratio2) //PMT 39mm
619 {
620 forb=1;
621 pathL=(m_scinLength/2+emtPos.z())/costheta;
622 }
623 }
624 else
625 {
626 if (G4UniformRand()<ratio3)
627 {
628 forb=0;
629 pathL=(1.5*m_scinLength+emtPos.z())/costheta;
630 }
631 }
632 }
633 else //Air
634 {
635 G4double tempran=G4UniformRand();
636 if (tempran<ratio3)
637 {
638 forb=0;
639 pathL=(1.5*m_scinLength+emtPos.z())/costheta;
640 }
641 else if (tempran>ratio1 && G4UniformRand()<ratio3)
642 {
643 forb=1;
644 pathL=(2.5*m_scinLength+emtPos.z())/costheta;
645 }
646 }
647 }
648 }
649
650 G4double convFactor=180./3.1415926;
651 if (theta>asin(1)-thetaR)
652 {
653 G4double alpha = G4UniformRand()*asin(1.);
654 G4int nRef1 = pathL*sintheta*cos(alpha)/50.0+0.5;
655 G4int nRef2 = pathL*sintheta*sin(alpha)/58.9+0.5;
656 G4double beta1=acos(sintheta*cos(alpha));
657 G4double beta2=acos(sintheta*sin(alpha));
658 beta2 -= 3.75/convFactor;
659 G4double R21,R22;
660 R21=Reflectivity(m_refIndex,1.0,beta1);
661 R22=Reflectivity(m_refIndex,1.0,beta2);
662 for (G4int i=0;i<nRef1;i++)
663 {
664 if (G4UniformRand()<(1-R21) && G4UniformRand()<0.15)
665 pathL=-9;
666 }
667 for (G4int i=0;i<nRef2;i++)
668 {
669 if (G4UniformRand()<(1-R22) && G4UniformRand()<0.15)
670 pathL=-9;
671 }
672 }
673}
674
675
676//void BesTofDigitizerBrV2::DirectPh(G4int partId, G4ThreeVector emtPos, G4double& pathL, G4int& forb)
677//{
678// //generation photon have random direction
679// //optical parameters of scintillation simulation
680// G4double cos_span=1-cos( asin(1./m_refIndex) );
681// G4double dcos, ran;
682// ran=G4UniformRand();
683// //assuming uniform phi distribution, simulate cos distr only
684// dcos=cos_span*(ran*2.0-1.0);
685// if (dcos > 0.0&& dcos != 1)
686// {
687// forb=0; //forward
688// pathL=(m_scinLength/2-emtPos.z())/(1.0-dcos);
689// }
690// else if (dcos < 0 && dcos != -1)
691// {
692// forb=1;
693// pathL=(m_scinLength/2+emtPos.z())/(1.0+dcos);
694// }
695//}
696
698{
699 G4double tmp_tauRatio,tmp_tau1,tmp_tau2,tmp_tau3;
700 tmp_tauRatio = m_tauRatio;
701 tmp_tau1 = m_tau1;
702 tmp_tau2 = m_tau2;
703 tmp_tau3 = m_tau3;
704 G4double UniformR = tmp_tauRatio/(1+tmp_tauRatio);
705 G4double EmissionTime;
706 if (G4UniformRand()>UniformR) {
707 while (1) {
708 EmissionTime = -tmp_tau2*log( G4UniformRand() );
709 if (G4UniformRand()-exp(EmissionTime/tmp_tau2-EmissionTime/tmp_tau1)>1.E-8)
710 break;
711 }
712 }
713 else EmissionTime = -tmp_tau3*log( G4UniformRand() );
714 return EmissionTime;
715}
716
717/*G4double BesTofDigitizerBrV2::GenPhoton(G4int partId)
718 {
719//get scintilation time
720G4double scinTime;
721static G4double scinA[26000];
722G4double random[2];
723G4double inverseRegion, ran1, ran2, problive;
724static G4int istore=-1;
725if(istore<0)
726{
727istore=1;
728for(G4int i=0;i<26000;i++)
729scinA[i]=Scintillation( (i+1) *0.001 , partId);
730}
731while(1)
732{
733HepRandom::getTheEngine()->flatArray(2, random);
734inverseRegion=random[0]*2.00975218688231;
735ran1=-4.5*log(1-inverseRegion/(0.45*4.5));
736ran2=0.45*exp(-1.*ran1/4.5)*random[1];
737problive=scinA[G4int(ran1*1000)+1];
738if(problive>=ran2)
739{ scinTime=ran1; break; }
740}
741return scinTime;
742}*/
743
745{
746 //get time transit spread
747 return G4RandGauss::shoot(m_ttsMean,m_ttsSigma);
748}
749
750void BesTofDigitizerBrV2::AccuSignal(G4double endTime, G4int forb)
751{
752 G4int ihst;
753 ihst=G4int(endTime/m_timeBinSize);
754 if (ihst>0 &&ihst<m_profBinN)
755 {
756 m_nPhot[ihst][forb]=m_nPhot[ihst][forb]+1;
757 m_totalPhot[forb]=m_totalPhot[forb]+1;
758 }
759}
760
762{
763 //std::cout << "BesTofDigitizerBrV2::TofPmtRspns()" << std::endl;
765 //to generate PMT signal shape for single photoelectron.
766 //only one time for a job.
767 G4double snpe[m_snpeBinN][2];
768
769 //Model: f(t)=Gain*mv_1pe* t**2 * exp-(t/tau)**2/normal-const
770 //normalization const =sqrt(pi)*tau*tau*tau/4.0
771 //G4double tau = m_riseTime;
772 G4double norma_const;
773 G4double echarge=1.6e-7; //in unit of PC
774
775 //time profile of pmt signals for Back and Forward PMT.
776 G4double profPmt[m_profBinN][2];
777
778 G4double t;
779 G4double tau;
780 G4int n1, n2, ii;
781 G4int phtn;
782 // forb = 0, east
783 for (G4int i=0;i<m_snpeBinN;i++)
784 {
785 tau = tofPara->GetBrERiseTime(scinNb);
786 norma_const = sqrt(CLHEP::pi)*tau*tau*tau/4.0; //in unit of ns**3
787 t = (i+1)*m_timeBinSize;
788 snpe[i][0] = m_Ce*t*t*exp(- (t/tau) * (t/tau) )/norma_const;// Pulse of one single photoelectron
789 }
790 // forb = 1, west
791 for (G4int i=0;i<m_snpeBinN;i++)
792 {
793 tau = tofPara->GetBrWRiseTime(scinNb);
794 norma_const = sqrt(CLHEP::pi)*tau*tau*tau/4.0; //in unit of ns**3
795 t = (i+1)*m_timeBinSize;
796 snpe[i][1] = m_Ce*t*t*exp(- (t/tau) * (t/tau) )/norma_const;
797 }
798 //for barrel and endcap tof: fb=2 or 1
799 G4int fb=2;
800 G4double Npoisson;
801 G4double pmtGain0,pmtGain,relativeGain,smearPMTgain;
802 G4double smearADC[2] = {0.,0.};
803 G4int runId = m_RealizationSvc->getRunId();
804 pmtGain0 = m_tofSimSvc->BarPMTGain();
805
806// if(runId>=-80000 && runId<=-9484)
807// {
808// // After TOF HV adjustment, PMT Gain of barrel TOF in MC is set to 5.E5
809// // High/Low Threshold for barrel: 500/125 mV
810// pmtGain0 = 6.E5;
811// }
812// else
813// {
814// // Before TOF HV adjustment, PMT Gain of barrel TOF in MC is set to 2.5E5
815// // High/Low Threshold for barrel: 600/150 mV
816// pmtGain0 = 5.E5;
817// }
818
819 G4double timeSmear = G4RandGauss::shoot(0,0.020);
820 if (runId>=-22913 && runId<=-20448) {//for 2011-psipp(20448-22913), smear barrel TOF resolution to ~78ps
821 timeSmear = G4RandGauss::shoot(0,0.040);
822 }
823 else if( ( runId>=-11396 && runId<=-8093 ) || ( runId>-80000 && runId<=-23463 ) ) {
824 timeSmear = G4RandGauss::shoot(0,0.025);
825 }
826
827 for (G4int j=0; j<fb; j++)
828 {
829 if (m_totalPhot[j] > 0)
830 {
831 n1=G4int(m_t1st[j]/m_timeBinSize);
832 n2=G4int(m_tLast[j]/m_timeBinSize);
833 //std::cout << "n1: " << n1 << std::endl;
834
835 for (G4int i=0;i<m_profBinN;i++)
836 profPmt[i][j]=0.0;
837
838 //generate PMT pulse
840 for (G4int i=n1;i<n2;i++)
841 {
842 phtn=m_nPhot[i][j];
843 if (phtn>0)
844 {
845 while(1) {
846 Npoisson = G4Poisson(10.0);
847 if(Npoisson>0.) break;
848 }
849 while(1) {
850 //pmtGain = j ? tofPara->GetBrWPMTgain(scinNb) : tofPara->GetBrEPMTgain(scinNb);
851 //relativeGain = j ? m_tofCaliSvc->BGainBackward(scinNb) : m_tofCaliSvc->BGainForward(scinNb);
852 relativeGain = j ? m_tofSimSvc->BarGain2(scinNb) : m_tofSimSvc->BarGain1(scinNb);
853 pmtGain = pmtGain0*relativeGain;
854 smearPMTgain = G4RandGauss::shoot(pmtGain,pmtGain/sqrt(Npoisson));
855 //smearPMTgain = pmtGain;
856 if(smearPMTgain>0) break;
857 }
858 smearADC[j] += phtn*smearPMTgain;
859
860 for (G4int ihst=0; ihst<m_snpeBinN; ihst++)
861 {
862 ii=i+ihst;
863 if (ii<m_profBinN)
864 profPmt[ii][j] += smearPMTgain*phtn*snpe[ihst][j];
865 else
866 break;
867 }
868 }
869 }
870
871 //add preamplifier and noise
872 for (int i=0;i<m_profBinN;i++)
873 {
874 if (profPmt[i][j]>0)
875 profPmt[i][j] = m_preGain*profPmt[i][j]+G4RandGauss::shoot(0,m_noiseSigma);
876 }
877
878 //get pulse height
879 G4double max=0;
880 for (int i=n1;i<m_profBinN;i++)
881 {
882 if (profPmt[i][j]>max)
883 max=profPmt[i][j];
884 }
885 if (m_G4Svc->TofRootFlag())
886 {
887 if (j==0) m_max0=max;
888 else m_max1=max;
889 }
890
891 G4double tmp_HLthresh, tmp_LLthresh, adcFactor;
892 G4double ratio;
893
894 if(runId>=-80000 && runId<=-9484)
895 {
896 // After TOF HV adjustment, PMT Gain of barrel TOF in MC is set to 5.E5
897 // High/Low Threshold for barrel: 500/125 mV
898 tmp_HLthresh = m_HLthresh;
899 tmp_LLthresh = m_LLthresh;
900 adcFactor = 5.89;
901 }
902 else
903 {
904 // Before TOF HV adjustment, PMT Gain of barrel TOF in MC is set to 2.5E5
905 // High/Low Threshold for barrel: 600/150 mV
906 tmp_HLthresh = 600;
907 tmp_LLthresh = 150;
908 adcFactor = 4.8;
909 }
910// if(runId>=-80000 && runId<=-9484)
911// {
912// ratio=2.16*1923.8/2197.8;
913// }
914// else
915// {
916// ratio = 2.11*2437.0/3102.9;
917// }
918 tmp_HLthresh = m_tofSimSvc->BarHighThres();
919 tmp_LLthresh = m_tofSimSvc->BarLowThres();
920 ratio = m_tofSimSvc->BarConstant();
921
922 //get final tdc and adc
923 if (max>=tmp_HLthresh)
924 {
925 for (int i=0;i<m_profBinN;i++)
926 {
927 if ( profPmt[i][j] >= tmp_LLthresh )
928 {
929 m_TDC[j] = i*m_timeBinSize + timeSmear +G4RandGauss::shoot(0,0.025); // Adding Electronical Uncertainty of 25ps
930// hhliu 20140724, setting tdc-smear for inner and outer separately
931 if( ( runId>=-11396 && runId<=-8093 ) || ( runId>-80000 && runId<=-23463 ) ) {
932 if( scinNb<88 ) {
933 m_TDC[j] = i*m_timeBinSize + timeSmear +G4RandGauss::shoot(0,0.055);
934 }
935 else {
936 m_TDC[j] = i*m_timeBinSize + timeSmear +G4RandGauss::shoot(0,0.062);
937 }
938 }
939
940 if( m_G4Svc->TofSaturationFlag())
941 {
942 //get ADC[j] using tofQElecSvc
943 double x = m_preGain*smearADC[j]*echarge*ratio;
944 if (j==0)
945 {
946 m_ADC[j] = m_tofQElecSvc->BQChannel1(scinNb,x);
947 }
948 else
949 {
950 m_ADC[j] = m_tofQElecSvc->BQChannel2(scinNb,x);
951 }
952 }
953 else
954 m_ADC[j] = m_preGain*smearADC[j]*echarge*adcFactor;
955
956
957 if (m_G4Svc->TofRootFlag())
958 {
959 if (j==0) {
960 m_tdc0 = m_TDC[0];
961 m_adc0 = m_ADC[0];
962 }
963 else {
964 m_tdc1 = m_TDC[1];
965 m_adc1 = m_ADC[1];
966 }
967 }
968 break;
969 }
970 }
971 }
972 }
973 }
974}
Double_t x[10]
EvtComplex exp(const EvtComplex &c)
Definition: EvtComplex.hh:252
const double alpha
double sin(const BesAngle a)
double cos(const BesAngle a)
G4TDigiCollection< BesTofDigi > BesTofDigitsCollection
G4THitsCollection< BesTofHit > BesTofHitsCollection
int n2
Definition: SD0Tag.cxx:55
int n1
Definition: SD0Tag.cxx:54
G4double Reflectivity(G4double n1, G4double n2, G4double theta)
void AccuSignal(G4double, G4int)
virtual void Digitize(ScintSingle *, BesTofDigitsCollection *)
G4double Scintillation(G4int)
void DirectPh(G4int, G4ThreeVector, G4double &, G4int &)
G4double BirksLaw(BesTofHit *hit)
void TofPmtAccum(BesTofHit *, G4int)
static BesTofGeoParameter * GetInstance()
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 BarAttenLength(unsigned int id)=0
virtual const double BarHighThres()=0
virtual const double BarGain1(unsigned int id)=0
int t()
Definition: t.c:1
#define ns(x)
Definition: xmltok.c:1504