2 /**************************************************************************
3 * Copyright(c) 1998-2000, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
8 * Permission to use, copy, modify and distribute this software and its *
9 * documentation strictly for non-commercial purposes is hereby granted *
10 * without fee, provided that the above copyright notice appears in all *
11 * copies and that both the copyright notice and this permission notice *
12 * appear in the supporting documentation. The authors make no claims *
13 * about the suitability of this software for any purpose. It is *
14 * provided "as is" without express or implied warranty. *
15 **************************************************************************/
19 /******************************************************************
20 * Produde digits from hits
21 * digits is TObject and includes
22 * We are writing array if C & A TDC
23 * C & A ADC (will need for slow simulation)
24 * TOF first particle C & A
25 * mean time and time difference (vertex position)
27 * Alla.Maevskaya@cern.ch
28 ****************************************************************/
40 #include "AliT0Digitizer.h"
43 #include "AliT0digit.h"
44 #include "AliRunDigitizer.h"
46 #include <AliLoader.h>
47 #include <AliRunLoader.h>
49 #include "AliT0Parameters.h"
51 ClassImp(AliT0Digitizer)
53 //___________________________________________
54 AliT0Digitizer::AliT0Digitizer() :AliDigitizer(),
58 ftimeCFD(new TArrayI(24)),
59 ftimeLED (new TArrayI(24)),
60 fADC(new TArrayI(24)),
61 fADC0 (new TArrayI(24)),
67 // Default ctor - don't use it
71 //___________________________________________
72 AliT0Digitizer::AliT0Digitizer(AliRunDigitizer* manager)
73 :AliDigitizer(manager),
77 ftimeCFD(new TArrayI(24)),
78 ftimeLED (new TArrayI(24)),
79 fADC(new TArrayI(24)),
80 fADC0 (new TArrayI(24)),
86 // ctor which should be used
88 AliDebug(1,"processed");
89 fParam = AliT0Parameters::Instance();
92 for (Int_t i=0; i<24; i++){
93 TGraph* gr = fParam ->GetAmpLED(i);
95 Int_t np = gr->GetN();
96 Double_t *x = gr->GetX();
97 Double_t *y = gr->GetY();
99 Double_t *x1 = new Double_t[np];
100 Double_t *y1 = new Double_t[np];
101 for (Int_t ii=0; ii<np; ii++) {
106 TGraph *grInverse = new TGraph(np,y1,x1);
107 fAmpLED.AddAtAndExpand(grInverse,i);
110 for (Int_t i=0; i<24; i++){
111 TGraph* grq = fParam ->GetQTC(i);
113 Int_t npq = grq->GetN();
114 Double_t *xq = grq->GetX();
115 Double_t *yq = grq->GetY();
116 Double_t *x1q = new Double_t[npq];
117 Double_t *y1q = new Double_t[npq];
118 for (Int_t ii=1; ii<npq; ii++) {
122 TGraph *grInverseQTC = new TGraph(npq,y1q,x1q);
123 fAmpQTC.AddAtAndExpand(grInverseQTC,i);
129 //------------------------------------------------------------------------
130 AliT0Digitizer::~AliT0Digitizer()
143 //------------------------------------------------------------------------
144 Bool_t AliT0Digitizer::Init()
151 //---------------------------------------------------------------------
152 void AliT0Digitizer::Exec(Option_t* /*option*/)
156 Produde digits from hits
157 digits is TObject and includes
158 We are writing array if C & A TDC
159 C & A ADC (will need for slow simulation)
160 TOF first particle C & A
161 mean time and time difference (vertex position)
166 AliRunLoader *outRL = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
167 AliLoader * pOutStartLoader = outRL->GetLoader("T0Loader");
169 AliDebug(1,"start...");
172 // From hits to digits
176 Int_t volume, pmt, trCFD, trLED;
181 Float_t time[24], besttime[24], timeGaus[24] ;
182 //Q->T-> coefficients !!!! should be asked!!!
183 Float_t timeDelayCFD[24];
184 Int_t threshold =50; //photoelectrons
185 Float_t zdetA, zdetC;
186 Int_t sumMultCoeff = 100;
189 Int_t ph2Mip = fParam->GetPh2Mip();
190 Float_t channelWidth = fParam->GetChannelWidth() ;
191 Float_t delayVertex = fParam->GetTimeDelayTVD();
194 zdetC = TMath::Abs(fParam->GetZPosition("T0/C/PMT1"));
195 zdetA = TMath::Abs(fParam->GetZPosition("T0/A/PMT15"));
197 // printf(" !!!!!Z det A = %f C = % f",zdetA,zdetC);
201 Int_t nFiles=fManager->GetNinputs();
202 for (Int_t inputFile=0; inputFile<nFiles; inputFile++) {
203 if (inputFile < nFiles-1) {
204 AliWarning(Form("ignoring input stream %d", inputFile));
209 Float_t besttimeC=99999.;
210 Float_t besttimeA=99999.;
211 Int_t pmtBestC=99999;
212 Int_t pmtBestA=99999;
213 Int_t timeDiff=99999, meanTime=99999;
214 Int_t sumMult =0; fSumMult=0;
215 bestATDC = 99999; bestCTDC = 99999;
222 for (Int_t i0=0; i0<24; i0++)
224 time[i0]=besttime[i0]=timeGaus[i0]=999999; countE[i0]=0;
226 AliRunLoader * inRL = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
227 AliLoader * pInStartLoader = inRL->GetLoader("T0Loader");
228 if (!inRL->GetAliRun()) inRL->LoadgAlice();
229 fT0 = (AliT0*)inRL ->GetAliRun()->GetDetector("T0");
232 pInStartLoader->LoadHits("READ");//probably it is necessary to load them before
233 fHits = fT0->Hits ();
234 TTree *th = pInStartLoader->TreeH();
235 brHits = th->GetBranch("T0");
237 fT0->SetHitsAddressBranch(brHits);
239 AliWarning("Branch T0 hit not found for this event");
240 // fT0->AddDigit(bestATDC,bestCTDC,meanTime,timeDiff,fSumMult, refpoint,
241 // ftimeCFD,fADC0,ftimeLED,fADC);
244 Int_t ntracks = (Int_t) th->GetEntries();
246 if (ntracks<=0) return;
247 // Start loop on tracks in the hits containers
248 for (Int_t track=0; track<ntracks;track++) {
249 brHits->GetEntry(track);
250 nhits = fHits->GetEntriesFast();
251 for (hit=0;hit<nhits;hit++)
253 startHit = (AliT0hit*) fHits->UncheckedAt(hit);
255 AliError("The unchecked hit doesn't exist");
260 Double_t e=startHit->Etot();
261 volume = startHit->Volume();
265 besttime[numpmt] = startHit->Time();
266 if(besttime[numpmt]<time[numpmt])
268 time[numpmt]=besttime[numpmt];
270 } //photoelectron accept
274 //spread time A&C by 25ps && besttime
275 Float_t c = 0.0299792; // cm/ps
277 Float_t koef=(zdetA-zdetC)/c; //correction position difference by cable
278 for (Int_t ipmt=0; ipmt<12; ipmt++){
279 if(countE[ipmt] > threshold) {
280 timeGaus[ipmt]=gRandom->Gaus(time[ipmt],25)+koef;
281 if(timeGaus[ipmt]<besttimeC){
282 besttimeC=timeGaus[ipmt]; //timeC
286 for ( Int_t ipmt=12; ipmt<24; ipmt++){
287 if(countE[ipmt] > threshold) {
288 timeGaus[ipmt]=gRandom->Gaus(time[ipmt],25);
289 if(timeGaus[ipmt]<besttimeA) {
290 besttimeA=timeGaus[ipmt]; //timeA
295 timeDelayCFD[0] = fParam->GetTimeDelayCFD(0);
297 for (Int_t i=0; i<24; i++)
299 Float_t al = countE[i];
300 if (al>threshold && timeGaus[i]<50000 ) {
303 // phe -> mV 0.3; 1MIP ->500phe -> ln (amp (mV)) = 5;
304 // max 200ns, HIJING mean 50000phe -> 15000mv -> ln = 15 (s zapasom)
306 qt= al/ph2Mip; // 50mv/Mip amp in mV
307 // before will we have calibration for high multiplicity
308 // if (qt > 115.) qt =115.; //must be fix!!!
310 timeDelayCFD[i] = fParam->GetTimeDelayCFD(i);
311 trCFD = Int_t (timeGaus[i]/channelWidth + timeDelayCFD[i]);
313 TGraph* gr = ((TGraph*)fAmpLED.At(i));
314 if(gr) sl = gr->Eval(qt);
316 TGraph* gr1 = ((TGraph*)fAmpQTC.At(i));
317 if(gr1) qtCh = gr1->Eval(qt);
320 fADC->AddAt(Int_t(qtCh),i);
321 // sumMult += Int_t ((al*gain[i]/ph2Mip)*50) ;
322 sumMult += Int_t (qtCh/sumMultCoeff) ;
325 TGraph *fu=(TGraph*) fParam ->GetWalk(i);
326 if(fu) slew=fu->Eval(Float_t(qtCh));
328 // trCFD=trCFD-Int_t(fMaxValue[i]-slew);
329 trCFD = trCFD + Int_t(slew); //for the same channel as cosmic
330 ftimeCFD->AddAt(Int_t (trCFD),i);
331 trLED = Int_t(trCFD + sl );
332 ftimeLED->AddAt(trLED,i);
333 AliDebug(1,Form(" pmt %i : delay %f time in ns %f time in channels %i LEd %i ", i, timeDelayCFD[i], timeGaus[i],trCFD, trLED ));
334 AliDebug(1,Form(" qt in MIP %f led-cfd in %f qt in channels %f ",qt, sl, qtCh));
339 //folding with alignmentz position distribution
340 if( besttimeC > 10000. && besttimeC <15000)
341 bestCTDC=Int_t ((besttimeC+timeDelayCFD[pmtBestC])
344 if( besttimeA > 10000. && besttimeA <15000)
345 bestATDC=Int_t ((besttimeA+timeDelayCFD[pmtBestA])
348 if (bestATDC < 99999 && bestCTDC < 99999)
350 timeDiff=Int_t (((besttimeC-besttimeA)+1000*delayVertex)
352 meanTime=Int_t (((besttimeC+besttimeA)/2. )/channelWidth);
355 if (sumMult > threshold){
356 fSumMult = Int_t (1000.* TMath::Log(Double_t(sumMult) / Double_t(sumMultCoeff))
358 AliDebug(10,Form("summult mv %i mult in chammens %i in ps %f ",
359 sumMult, fSumMult, fSumMult*channelWidth));
362 fT0->AddDigit(bestATDC,bestCTDC,meanTime,timeDiff,fSumMult, refpoint,
363 ftimeCFD,fADC0,ftimeLED,fADC);
366 AliDebug(10,Form(" Digits wrote refpoint %i bestATDC %i bestCTDC %i meanTime %i timeDiff %i fSumMult %i ",refpoint ,bestATDC,bestCTDC,meanTime,timeDiff,fSumMult ));
367 pOutStartLoader->UnloadHits();
368 } //input streams loop
371 pOutStartLoader->LoadDigits("UPDATE");
372 TTree *treeD = pOutStartLoader->TreeD();
374 pOutStartLoader->MakeTree("D");
375 treeD = pOutStartLoader->TreeD();
378 fT0 = (AliT0*)outRL ->GetAliRun()->GetDetector("T0");
379 // Make a branch in the tree
380 fT0->MakeBranch("D");
383 pOutStartLoader->WriteDigits("OVERWRITE");
386 pOutStartLoader->UnloadDigits();