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 "AliDigitizationInput.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(AliDigitizationInput* digInput)
73 :AliDigitizer(digInput),
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++) {
105 TGraph *grInverse = new TGraph(np,y1,x1);
106 fAmpLED.AddAtAndExpand(grInverse,i);
107 if (x1) delete [] x1;
108 if (y1) delete [] y1;
111 for (Int_t i=0; i<24; i++){
112 TGraph* grq = fParam ->GetQTC(i);
114 Int_t npq = grq->GetN();
115 Double_t *xq = grq->GetX();
116 Double_t *yq = grq->GetY();
117 Double_t *x1q = new Double_t[npq];
118 Double_t *y1q = new Double_t[npq];
119 for (Int_t ii=1; ii<npq; ii++) {
123 TGraph *grInverseQTC = new TGraph(npq,y1q,x1q);
124 fAmpQTC.AddAtAndExpand(grInverseQTC,i);
125 if (x1q) delete [] x1q;
126 if (y1q) delete [] y1q;
133 //------------------------------------------------------------------------
134 AliT0Digitizer::~AliT0Digitizer()
147 //------------------------------------------------------------------------
148 Bool_t AliT0Digitizer::Init()
155 //---------------------------------------------------------------------
156 void AliT0Digitizer::Digitize(Option_t* /*option*/)
160 Produde digits from hits
161 digits is TObject and includes
162 We are writing array if C & A TDC
163 C & A ADC (will need for slow simulation)
164 TOF first particle C & A
165 mean time and time difference (vertex position)
170 AliRunLoader *outRL = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName());
171 AliLoader * pOutStartLoader = outRL->GetLoader("T0Loader");
173 AliDebug(1,"start...");
176 // From hits to digits
180 Int_t volume, pmt, trCFD, trLED;
185 Float_t time[24], besttime[24], timeGaus[24] ;
186 //Q->T-> coefficients !!!! should be asked!!!
187 Float_t timeDelayCFD[24];
188 Int_t threshold =50; //photoelectrons
189 Float_t zdetA, zdetC;
190 Int_t sumMultCoeff = 100;
193 Int_t ph2Mip = fParam->GetPh2Mip();
194 Float_t channelWidth = fParam->GetChannelWidth() ;
195 Float_t delayVertex = fParam->GetTimeDelayTVD();
198 zdetC = TMath::Abs(fParam->GetZPosition("T0/C/PMT1"));
199 zdetA = TMath::Abs(fParam->GetZPosition("T0/A/PMT15"));
201 // printf(" !!!!!Z det A = %f C = % f",zdetA,zdetC);
205 Int_t nFiles=fDigInput->GetNinputs();
206 for (Int_t inputFile=0; inputFile<nFiles; inputFile++) {
207 if (inputFile < nFiles-1) {
208 AliWarning(Form("ignoring input stream %d", inputFile));
213 Float_t besttimeC=99999.;
214 Float_t besttimeA=99999.;
215 Int_t pmtBestC=99999;
216 Int_t pmtBestA=99999;
217 Int_t timeDiff=99999, meanTime=99999;
218 Int_t sumMult =0; fSumMult=0;
219 bestATDC = 99999; bestCTDC = 99999;
226 for (Int_t i0=0; i0<24; i0++)
228 time[i0]=besttime[i0]=timeGaus[i0]=999999; countE[i0]=0;
230 AliRunLoader * inRL = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(inputFile));
231 AliLoader * pInStartLoader = inRL->GetLoader("T0Loader");
232 if (!inRL->GetAliRun()) inRL->LoadgAlice();
233 fT0 = (AliT0*)inRL ->GetAliRun()->GetDetector("T0");
236 pInStartLoader->LoadHits("READ");//probably it is necessary to load them before
237 fHits = fT0->Hits ();
238 TTree *th = pInStartLoader->TreeH();
239 brHits = th->GetBranch("T0");
241 fT0->SetHitsAddressBranch(brHits);
243 AliWarning("Branch T0 hit not found for this event");
244 // fT0->AddDigit(bestATDC,bestCTDC,meanTime,timeDiff,fSumMult, refpoint,
245 // ftimeCFD,fADC0,ftimeLED,fADC);
248 Int_t ntracks = (Int_t) th->GetEntries();
250 if (ntracks<=0) return;
251 // Start loop on tracks in the hits containers
252 for (Int_t track=0; track<ntracks;track++) {
253 brHits->GetEntry(track);
254 nhits = fHits->GetEntriesFast();
255 for (hit=0;hit<nhits;hit++)
257 startHit = (AliT0hit*) fHits->UncheckedAt(hit);
259 AliError("The unchecked hit doesn't exist");
264 Double_t e=startHit->Etot();
265 volume = startHit->Volume();
269 besttime[numpmt] = startHit->Time();
270 if(besttime[numpmt]<time[numpmt])
272 time[numpmt]=besttime[numpmt];
274 } //photoelectron accept
278 //spread time A&C by 25ps && besttime
279 Float_t c = 0.0299792; // cm/ps
281 Float_t koef=(zdetA-zdetC)/c; //correction position difference by cable
282 for (Int_t ipmt=0; ipmt<12; ipmt++){
283 if(countE[ipmt] > threshold) {
284 timeGaus[ipmt]=gRandom->Gaus(time[ipmt],25)+koef;
285 if(timeGaus[ipmt]<besttimeC){
286 besttimeC=timeGaus[ipmt]; //timeC
290 for ( Int_t ipmt=12; ipmt<24; ipmt++){
291 if(countE[ipmt] > threshold) {
292 timeGaus[ipmt]=gRandom->Gaus(time[ipmt],25);
293 if(timeGaus[ipmt]<besttimeA) {
294 besttimeA=timeGaus[ipmt]; //timeA
299 timeDelayCFD[0] = fParam->GetTimeDelayCFD(0);
301 for (Int_t i=0; i<24; i++)
303 Float_t al = countE[i];
304 if (al>threshold && timeGaus[i]<50000 ) {
307 // phe -> mV 0.3; 1MIP ->500phe -> ln (amp (mV)) = 5;
308 // max 200ns, HIJING mean 50000phe -> 15000mv -> ln = 15 (s zapasom)
310 qt= al/ph2Mip; // 50mv/Mip amp in mV
311 // before will we have calibration for high multiplicity
312 // if (qt > 115.) qt =115.; //must be fix!!!
314 timeDelayCFD[i] = fParam->GetTimeDelayCFD(i);
315 trCFD = Int_t (timeGaus[i]/channelWidth + timeDelayCFD[i]);
317 TGraph* gr = ((TGraph*)fAmpLED.At(i));
318 if(gr) sl = gr->Eval(qt);
320 TGraph* gr1 = ((TGraph*)fAmpQTC.At(i));
321 if(gr1) qtCh = gr1->Eval(qt);
324 fADC->AddAt(Int_t(qtCh),i);
325 // sumMult += Int_t ((al*gain[i]/ph2Mip)*50) ;
326 sumMult += Int_t (qtCh/sumMultCoeff) ;
329 TGraph *fu=(TGraph*) fParam ->GetWalk(i);
330 if(fu) slew=fu->Eval(Float_t(qtCh));
332 // trCFD=trCFD-Int_t(fMaxValue[i]-slew);
333 trCFD = trCFD + Int_t(slew); //for the same channel as cosmic
334 ftimeCFD->AddAt(Int_t (trCFD),i);
335 trLED = Int_t(trCFD + sl );
336 ftimeLED->AddAt(trLED,i);
337 AliDebug(1,Form(" pmt %i : delay %f time in ns %f time in channels %i LEd %i ", i, timeDelayCFD[i], timeGaus[i],trCFD, trLED ));
338 AliDebug(1,Form(" qt in MIP %f led-cfd in %f qt in channels %f ",qt, sl, qtCh));
343 //folding with alignmentz position distribution
344 if( besttimeC > 10000. && besttimeC <15000)
345 bestCTDC=Int_t ((besttimeC+timeDelayCFD[pmtBestC])
348 if( besttimeA > 10000. && besttimeA <15000)
349 bestATDC=Int_t ((besttimeA+timeDelayCFD[pmtBestA])
352 if (bestATDC < 99999 && bestCTDC < 99999)
354 timeDiff=Int_t (((besttimeC-besttimeA)+1000*delayVertex)
356 meanTime=Int_t (((besttimeC+besttimeA)/2. )/channelWidth);
359 if (sumMult > threshold){
360 fSumMult = Int_t (1000.* TMath::Log(Double_t(sumMult) / Double_t(sumMultCoeff))
362 AliDebug(10,Form("summult mv %i mult in chammens %i in ps %f ",
363 sumMult, fSumMult, fSumMult*channelWidth));
366 fT0->AddDigit(bestATDC,bestCTDC,meanTime,timeDiff,fSumMult, refpoint,
367 ftimeCFD,fADC0,ftimeLED,fADC);
370 AliDebug(10,Form(" Digits wrote refpoint %i bestATDC %i bestCTDC %i meanTime %i timeDiff %i fSumMult %i ",refpoint ,bestATDC,bestCTDC,meanTime,timeDiff,fSumMult ));
371 pOutStartLoader->UnloadHits();
372 } //input streams loop
375 pOutStartLoader->LoadDigits("UPDATE");
376 TTree *treeD = pOutStartLoader->TreeD();
378 pOutStartLoader->MakeTree("D");
379 treeD = pOutStartLoader->TreeD();
382 fT0 = (AliT0*)outRL ->GetAliRun()->GetDetector("T0");
383 // Make a branch in the tree
384 fT0->MakeBranch("D");
387 pOutStartLoader->WriteDigits("OVERWRITE");
390 pOutStartLoader->UnloadDigits();