00001
00008 #include "StBbcSimulationMaker.h"
00009 #include "g2t/St_g2t_bbc_Module.h"
00010 #include "TDataSetIter.h"
00011 #include "TRandom.h"
00012 #include <vector>
00013 #include "StEvent.h"
00014 #include "StTriggerDetectorCollection.h"
00015 #include "StBbcTriggerDetector.h"
00016 #include "StMessMgr.h"
00017
00018 TRandom BbcRndm = TRandom(0);
00019
00020 ClassImp(StBbcSimulationMaker)
00021
00022 const float BbcTimingRMS = 900.E-12;
00029 const u_short NPMTsmall1 = 16;
00030 const u_short NPMTlarge1 = 8;
00031 const u_short NPMT1 = NPMTsmall1+NPMTlarge1;
00032 const u_short NPMT2 = 2*NPMT1;
00033 const float dE1MIPper_gcm2 = 1.95E-3;
00034 const float PolystereneDensity = 1.032;
00035 const float TyleThickness = 1.;
00036 const float dE_1MIP = dE1MIPper_gcm2*PolystereneDensity*TyleThickness;
00037 const float NPhotoelectrons_1MIP = 15.;
00038 const float pC_per_Photoelectron = 0.3;
00039 const short NADCbins = 256;
00040 const short NTDCbins = 256;
00041 const float pC_perADCbin = 0.25;
00042
00043
00044 const float ADC0 = 0.;
00045 const float s_perTDCbin = .1E-9;
00046 const float TDC0 = 0.;
00047 const float OuterFactor = 0.8;
00048
00049
00050
00051 const float SinglePhotoElectronResolution = 0.3;
00052
00053
00054
00055
00056
00057 bool IsSmall(short iPMT)
00058 {
00060 if ( 0<=iPMT && iPMT<NPMTsmall1) return 1;
00061 if ( NPMT1<=iPMT && iPMT<NPMT1+NPMTsmall1) return 1;
00062 return 0;
00063 }
00064
00065
00066 class BbcTOF
00067 {
00072 private:
00073 vector<float> Times;
00074
00075 public:
00076 BbcTOF():Times(vector<float>(NPMT2)){};
00077 ~BbcTOF(){};
00078 void AddTOF(u_short ipmt, float time)
00079 {
00083 if (Times[ipmt]==0 || Times[ipmt]>time) {Times[ipmt]=time;}
00084 }
00085 float GetTOF(u_short ipmt)
00086 {
00088
00089 float tof = 0;
00090 if (Times[ipmt]!=0.){ tof = Times[ipmt]+BbcRndm.Gaus(0.,BbcTimingRMS); }
00091 return tof;
00092 }
00093 short GetTDC(u_short ipmt)
00094 {
00096 float T = this->GetTOF(ipmt);
00097
00098 if (T<TDC0) {return 0;}
00099 short N = (short)((T-TDC0)/s_perTDCbin);
00100 if (N>=NTDCbins) {return NTDCbins-1;}
00101 return N;
00102 }
00103
00104 };
00105
00106 class BbcDE
00107 {
00113 private:
00114 vector<float> dE;
00115
00116 public:
00117 BbcDE():dE(vector<float>(NPMT2)){};
00118 ~BbcDE(){};
00119 void AddDE(u_short ipmt, float de)
00120 {
00121 if (!IsSmall(ipmt)) {de *= OuterFactor;}
00122 dE[ipmt] += de;
00123 }
00124 float GetDE(u_short ipmt)
00125 {
00127
00128 float PoissonMean = dE[ipmt]/dE_1MIP*NPhotoelectrons_1MIP;
00129 short NPhotoelectrons = BbcRndm.Poisson(PoissonMean);
00130 float Q = pC_per_Photoelectron*
00131 (1+BbcRndm.Gaus(0.,SinglePhotoElectronResolution))*
00132 NPhotoelectrons;
00133 return Q;
00134 }
00135 short GetADC(u_short ipmt)
00136 {
00138 float A = this->GetDE(ipmt);
00139 if (A<ADC0) {return 0;}
00140 short N = (short)((A-ADC0)/pC_perADCbin);
00141 if (N>=NADCbins) {return NADCbins-1;}
00142 return N;
00143 }
00144 };
00145
00147
00154 StBbcSimulationMaker::StBbcSimulationMaker(const char *name):StMaker(name)
00155 {
00158 }
00159
00160
00162
00171 StBbcSimulationMaker::~StBbcSimulationMaker()
00172 {
00173
00174 }
00175
00176
00177
00179 Int_t StBbcSimulationMaker::Init()
00180 {
00203
00204
00205 for (short iew=1; iew<=2; iew++)
00206 {
00207 short EW = iew*1000;
00208
00209 short EWshift = (2-iew)*NPMT1;
00216
00217 Geant2PMT[EW+113] = EWshift+1;
00218
00219 Geant2PMT[EW+123] = EWshift+2;
00220
00221 Geant2PMT[EW+133] = EWshift+3;
00222
00223 Geant2PMT[EW+143] = EWshift+4;
00224
00225 Geant2PMT[EW+153] = EWshift+5;
00226
00227 Geant2PMT[EW+163] = EWshift+6;
00228
00229 Geant2PMT[EW+111] = EWshift+7;
00230
00231 Geant2PMT[EW+112] = EWshift+8;
00232
00233 Geant2PMT[EW+121] = EWshift+7;
00234
00235 Geant2PMT[EW+122] = EWshift+9;
00236
00237 Geant2PMT[EW+131] = EWshift+10;
00238
00239 Geant2PMT[EW+132] = EWshift+11;
00240
00241 Geant2PMT[EW+141] = EWshift+12;
00242
00243 Geant2PMT[EW+142] = EWshift+13;
00244
00245 Geant2PMT[EW+151] = EWshift+12;
00246
00247 Geant2PMT[EW+152] = EWshift+14;
00248
00249 Geant2PMT[EW+161] = EWshift+15;
00250
00251 Geant2PMT[EW+162] = EWshift+16;
00252
00253
00254 Geant2PMT[EW+213] = EWshift+17;
00255
00256 Geant2PMT[EW+223] = EWshift+18;
00257
00258 Geant2PMT[EW+233] = EWshift+18;
00259
00260 Geant2PMT[EW+243] = EWshift+19;
00261
00262 Geant2PMT[EW+253] = EWshift+20;
00263
00264 Geant2PMT[EW+263] = EWshift+20;
00265
00266 Geant2PMT[EW+211] = EWshift+21;
00267
00268 Geant2PMT[EW+212] = EWshift+21;
00269
00270 Geant2PMT[EW+221] = EWshift+21;
00271
00272 Geant2PMT[EW+222] = EWshift+22;
00273
00274 Geant2PMT[EW+231] = EWshift+22;
00275
00276 Geant2PMT[EW+232] = EWshift+22;
00277
00278 Geant2PMT[EW+241] = EWshift+23;
00279
00280 Geant2PMT[EW+242] = EWshift+23;
00281
00282 Geant2PMT[EW+251] = EWshift+23;
00283
00284 Geant2PMT[EW+252] = EWshift+24;
00285
00286 Geant2PMT[EW+261] = EWshift+24;
00287
00288 Geant2PMT[EW+262] = EWshift+24;
00289 }
00290
00291
00292 #ifdef BbcSimQa
00293 QaFile = new TFile("StBbcSimQa.root","recreate");
00294 QaBbcPmtdE = new TH1F("QaBbcPmtdE","BBC PMT",NPMT2,-0.5,NPMT2-0.5);
00295 QaBbcPmtTime = new TH1F("QaBbcPmtTime","BBC PMT",NPMT2,-0.5,NPMT2-0.5);
00296 QaBbcEastVid =
00297 new TH1F("QaBbcEastVid","BBC PMT with East VID",256,-0.5,255.5);
00298 QaBbcWestVid =
00299 new TH1F("QaBbcWestVid","BBC PMT with West VID",256,-0.5,255.5);
00300 QaBbcEastPmt =
00301 new TH1F("QaBbcEastPmt","BBC PMT with East #", 256,-0.5,255.5);
00302 QaBbcWestPmt =
00303 new TH1F("QaBbcWestPmt","BBC PMT with West #", 256,-0.5,255.5);
00304
00305 typedef map<short,short>::const_iterator CI;
00306
00307 for (CI I=Geant2PMT.begin(); I!=Geant2PMT.end(); ++I)
00308 {
00309 PMT2Geant[(*I).second] = (*I).first;
00310 }
00311 #endif
00312
00313 return StMaker::Init();
00314 }
00315
00316
00317 Int_t StBbcSimulationMaker::Make()
00318 {
00320
00321 TDataSet* ds= GetInputDS("geant");
00322 assert(ds);
00323 StEvent* event = (StEvent*)GetInputDS("StEvent");
00324 assert(event);
00325 St_g2t_ctf_hit* g2t_bbc_hit = (St_g2t_ctf_hit*)ds->Find("g2t_bbc_hit");
00326
00327 if (g2t_bbc_hit)
00328 {
00329 short nBBChits = g2t_bbc_hit->GetNRows();
00330 BbcTOF TOFdata;
00331 BbcDE DEdata;
00332 g2t_ctf_hit_st *bbc_hit = g2t_bbc_hit->GetTable();
00333 for (short iBBChit=0; iBBChit<nBBChits; iBBChit++)
00334 {
00335 float De = bbc_hit[iBBChit].de;
00336 float TOF = bbc_hit[iBBChit].tof;
00337 short Vid = bbc_hit[iBBChit].volume_id;
00338
00339 short PMTid = Geant2PMT[Vid];
00340 if (PMTid == 0) {
00341 LOG_ERROR << "Cannot find a PMTid in Geant2PMT for Vid = " << Vid << endm;
00342 continue;
00343 }
00344 PMTid -= 1;
00345 DEdata.AddDE(PMTid,De);
00346 TOFdata.AddTOF(PMTid,TOF);
00347 }
00348
00349 StTriggerDetectorCollection* myTrig =event->triggerDetectorCollection();
00350 if (!myTrig) {
00351 Warning("Make"," NoStTriggerDetectorCollection, Make the new one\n");
00352 myTrig = new StTriggerDetectorCollection;
00353 event->setTriggerDetectorCollection(myTrig);
00354 }
00355 StBbcTriggerDetector& myBbc = myTrig->bbc();
00356
00357
00358 for (u_short iPMT = 0; iPMT<NPMT2; iPMT++)
00359 {
00360 short ADC = DEdata.GetADC(iPMT);
00361 #ifdef BbcSimQa
00362
00363 short Vid = PMT2Geant[iPMT+1];
00364
00365 if (Vid<2000) {QaBbcWestVid->Fill(ADC);}
00366 if (Vid>2000) {QaBbcEastVid->Fill(ADC);}
00367 if (iPMT<NPMT1) {QaBbcEastPmt->Fill(ADC);}
00368 if (NPMT1<=iPMT && iPMT<NPMT2) {QaBbcWestPmt->Fill(ADC);}
00369 QaBbcPmtTime->Fill(iPMT,TOFdata.GetTOF(iPMT));
00370 QaBbcPmtdE->Fill(iPMT,DEdata.GetDE(iPMT));
00371 #endif
00372 myBbc.setAdc(iPMT, ADC);
00373 if (IsSmall(iPMT))
00374 {
00375 short TDC = TOFdata.GetTDC(iPMT);
00376 myBbc.setTdc(iPMT, TDC);
00377 }
00378 }
00379 }
00380 else
00381 {
00382 gMessMgr->Info
00383 ("MLK StBbcSimulationMaker::Make() could not inst g2t_bbc_hit\n");
00384 }
00385
00386 return kStOK;
00387 }
00389 Int_t StBbcSimulationMaker::Finish()
00390 {
00392 #ifdef BbcSimQa
00393 QaFile->Write();
00394 QaFile->Close();
00395 #endif
00396 return StMaker::Finish();
00397 }
00398
00399
00400
00401
00402
00403
00404
00405
00406