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)),
68 // Default ctor - don't use it
72 //___________________________________________
73 AliT0Digitizer::AliT0Digitizer(AliRunDigitizer* manager)
74 :AliDigitizer(manager),
78 ftimeCFD(new TArrayI(24)),
79 ftimeLED (new TArrayI(24)),
80 fADC(new TArrayI(24)),
81 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 ->GetAmpLEDRec(i);
94 Int_t np = gr->GetN();
95 Double_t *x = gr->GetX();
96 Double_t *y = gr->GetY();
98 Double_t *x1 = new Double_t[np];
99 Double_t *y1 = new Double_t[np];
100 for (Int_t ii=0; ii<np; ii++) {
101 y1[ii]=y[np-ii-1]; x1[ii]=x[np-ii-1];
104 TGraph *grInverse = new TGraph(np,y1,x1);
107 fAmpLED.AddAtAndExpand(grInverse,i);
110 //------------------------------------------------------------------------
111 AliT0Digitizer::~AliT0Digitizer()
124 //------------------------------------------------------------------------
125 Bool_t AliT0Digitizer::Init()
132 //---------------------------------------------------------------------
133 void AliT0Digitizer::Exec(Option_t* /*option*/)
137 Produde digits from hits
138 digits is TObject and includes
139 We are writing array if C & A TDC
140 C & A ADC (will need for slow simulation)
141 TOF first particle C & A
142 mean time and time difference (vertex position)
147 AliRunLoader *outRL = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
148 AliLoader * pOutStartLoader = outRL->GetLoader("T0Loader");
150 AliDebug(1,"start...");
153 // From hits to digits
157 Int_t volume, pmt, trCFD, trLED;
159 Int_t bestATDC, bestCTDC, qtCh;
160 Float_t time[24], besttime[24], timeGaus[24] ;
161 //Q->T-> coefficients !!!! should be asked!!!
162 Float_t timeDelayCFD[24];
163 Int_t threshold =50; //photoelectrons
164 Float_t zdetA, zdetC;
165 Int_t sumMultCoeff = 100;
171 Int_t ph2Mip = fParam->GetPh2Mip();
172 Float_t channelWidth = fParam->GetChannelWidth() ;
173 Float_t delayVertex = fParam->GetTimeDelayTVD();
175 zdetC = TMath::Abs(fParam->GetZPosition("T0/C/PMT1"));
176 zdetA = TMath::Abs(fParam->GetZPosition("T0/A/PMT15"));
181 Int_t nFiles=fManager->GetNinputs();
182 for (Int_t inputFile=0; inputFile<nFiles; inputFile++) {
183 if (inputFile < nFiles-1) {
184 AliWarning(Form("ignoring input stream %d", inputFile));
189 Float_t besttimeC=99999.;
190 Float_t besttimeA=99999.;
193 Int_t timeDiff=999, meanTime=0;
194 Int_t sumMult =0; fSumMult=0;
195 bestATDC = 99999; bestCTDC = 99999;
202 for (Int_t i0=0; i0<24; i0++)
204 time[i0]=besttime[i0]=timeGaus[i0]=99999; countE[i0]=0;
206 AliRunLoader * inRL = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
207 AliLoader * pInStartLoader = inRL->GetLoader("T0Loader");
208 if (!inRL->GetAliRun()) inRL->LoadgAlice();
209 AliT0 *fT0 = (AliT0*)inRL ->GetAliRun()->GetDetector("T0");
212 pInStartLoader->LoadHits("READ");//probably it is necessary to load them before
213 TClonesArray *fHits = fT0->Hits ();
214 TTree *th = pInStartLoader->TreeH();
215 brHits = th->GetBranch("T0");
217 fT0->SetHitsAddressBranch(brHits);
219 AliError("Branch T0 hit not found");
222 Int_t ntracks = (Int_t) th->GetEntries();
224 if (ntracks<=0) return;
225 // Start loop on tracks in the hits containers
226 for (Int_t track=0; track<ntracks;track++) {
227 brHits->GetEntry(track);
228 nhits = fHits->GetEntriesFast();
229 for (hit=0;hit<nhits;hit++)
231 startHit = (AliT0hit*) fHits->UncheckedAt(hit);
233 AliError("The unchecked hit doesn't exist");
238 Double_t e=startHit->Etot();
239 volume = startHit->Volume();
241 // if(e>0 && RegisterPhotoE(numpmt,e)) {
244 besttime[numpmt] = startHit->Time();
245 if(besttime[numpmt]<time[numpmt])
247 time[numpmt]=besttime[numpmt];
249 } //photoelectron accept
253 //spread time A&C by 25ps && besttime
254 Float_t c = 0.0299792; // cm/ps
256 Float_t koef=(zdetA-zdetC)/c; //correction position difference by cable
257 for (Int_t ipmt=0; ipmt<12; ipmt++){
258 if(countE[ipmt] > threshold) {
259 timeGaus[ipmt]=gRandom->Gaus(time[ipmt],25)+koef;
260 if(timeGaus[ipmt]<besttimeC){
261 besttimeC=timeGaus[ipmt]; //timeC
265 for ( Int_t ipmt=12; ipmt<24; ipmt++){
266 if(countE[ipmt] > threshold) {
267 timeGaus[ipmt]=gRandom->Gaus(time[ipmt],25);
268 if(timeGaus[ipmt]<besttimeA) {
269 besttimeA=timeGaus[ipmt]; //timeA
274 timeDelayCFD[0] = fParam->GetTimeDelayCFD(0);
275 Int_t meanTimeDelay=200;
277 for (Int_t i=0; i<24; i++)
279 Float_t al = countE[i];
280 if (al>threshold && timeGaus[i]<50000 ) {
283 // phe -> mV 0.3; 1MIP ->500phe -> ln (amp (mV)) = 5;
284 // max 200ns, HIJING mean 50000phe -> 15000mv -> ln = 15 (s zapasom)
286 qt= 50.*al/ph2Mip; // 50mv/Mip amp in mV
288 timeDelayCFD[i] = fParam->GetTimeDelayCFD(i);
289 trCFD = Int_t (timeGaus[i]/channelWidth + (timeDelayCFD[i]-timeDelayCFD[0]));
290 TGraph* gr = ((TGraph*)fAmpLED.At(i));
293 trLED = Int_t(( timeGaus[i] + 1000*sl )/channelWidth);
294 qtCh=Int_t (1000.*TMath::Log(qt) / channelWidth);
297 ftimeLED->AddAt(trLED,i);
298 // sumMult += Int_t ((al*gain[i]/ph2Mip)*50) ;
299 sumMult += Int_t (qt/sumMultCoeff) ;
302 TGraph *fu=(TGraph*) fParam ->GetWalk(i) ;
303 Float_t slew=fu->Eval(Float_t(qtCh));
304 hr=fu->GetHistogram();
305 Float_t maxValue=hr->GetMaximum(50);
306 trCFD=trCFD-Int_t((maxValue-slew)/channelWidth);
307 ftimeCFD->AddAt(Int_t (trCFD),i);
308 AliDebug(10,Form(" pmt %i : time in ns %f time in channels %i ",
309 i, timeGaus[i],trCFD ));
310 AliDebug(10,Form(" qt in mV %f qt in ns %f qt in channels %i ",qt,
311 TMath::Log(qt), qtCh));
316 //folding with alignmentz position distribution
317 if( besttimeC > 10000. && besttimeC <15000)
318 bestCTDC=Int_t ((besttimeC+timeDelayCFD[pmtBestC])
321 if( besttimeA > 10000. && besttimeA <15000)
322 bestATDC=Int_t ((besttimeA+timeDelayCFD[pmtBestA])
325 if (bestATDC < 99999 && bestCTDC < 99999)
327 timeDiff=Int_t (((besttimeC-besttimeA)+1000*delayVertex)
329 meanTime=Int_t (((besttimeC+besttimeA)/2. + meanTimeDelay)/channelWidth);
331 AliDebug(10,Form(" time A& C %i %i time diff && mean time in channels %i %i",bestATDC,bestCTDC, timeDiff, meanTime));
333 if (sumMult > threshold){
334 fSumMult = Int_t (1000.* TMath::Log(Double_t(sumMult) / Double_t(sumMultCoeff))
336 AliDebug(10,Form("summult mv %i mult in chammens %i in ps %i ",
337 sumMult, fSumMult, fSumMult*channelWidth));
340 fT0->AddDigit(bestATDC,bestCTDC,meanTime,timeDiff,fSumMult, refpoint,
341 ftimeCFD,fADC0,ftimeLED,fADC);
344 AliDebug(10,Form(" Digits wrote refpoint %i bestATDC %i bestCTDC %i meanTime %i timeDiff %i fSumMult %i ",refpoint ,bestATDC,bestCTDC,meanTime,timeDiff,fSumMult ));
345 pOutStartLoader->UnloadHits();
346 } //input streams loop
349 pOutStartLoader->LoadDigits("UPDATE");
350 TTree *treeD = pOutStartLoader->TreeD();
352 pOutStartLoader->MakeTree("D");
353 treeD = pOutStartLoader->TreeD();
356 fT0 = (AliT0*)outRL ->GetAliRun()->GetDetector("T0");
357 // Make a branch in the tree
358 fT0->MakeBranch("D");
361 pOutStartLoader->WriteDigits("OVERWRITE");
364 pOutStartLoader->UnloadDigits();