]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0Digitizer.cxx
Updated list of classes
[u/mrichter/AliRoot.git] / T0 / AliT0Digitizer.cxx
1
2 /**************************************************************************
3  * Copyright(c) 1998-2000, ALICE Experiment at CERN, All rights reserved. *
4  *                                                                        *
5  * Author: The ALICE Off-line Project.                                    *
6  * Contributors are mentioned in the code where appropriate.              *
7  *                                                                        *
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  **************************************************************************/
16
17 /* $Id$ */
18
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)
26  *
27  *      Alla.Maevskaya@cern.ch 
28  ****************************************************************/
29
30
31 #include <TArrayI.h>
32 #include <TFile.h>
33 #include <TGraph.h>
34 #include <TH1F.h>
35 #include <TMath.h>
36 #include <TRandom.h>
37 #include <TTree.h> 
38
39 #include "AliLog.h"
40 #include "AliT0Digitizer.h"
41 #include "AliT0.h"
42 #include "AliT0hit.h"
43 #include "AliT0digit.h"
44 #include "AliRunDigitizer.h"
45 #include "AliRun.h"
46 #include <AliLoader.h>
47 #include <AliRunLoader.h>
48 #include <stdlib.h>
49 #include "AliT0Parameters.h"
50
51 ClassImp(AliT0Digitizer)
52
53 //___________________________________________
54   AliT0Digitizer::AliT0Digitizer()  :AliDigitizer(),
55                                      fT0(0),
56                                      fHits(0),
57                                      fdigits(0),
58                                      ftimeCFD(new TArrayI(24)), 
59                                      ftimeLED (new TArrayI(24)), 
60                                      fADC(new TArrayI(24)), 
61                                      fADC0 (new TArrayI(24)),
62                                      fSumMult(0),
63                                      fAmpLED(0),
64                                      fParam(0)
65
66
67 {
68 // Default ctor - don't use it
69   ;
70 }
71
72 //___________________________________________
73 AliT0Digitizer::AliT0Digitizer(AliRunDigitizer* manager) 
74   :AliDigitizer(manager),
75    fT0(0),
76    fHits(0),
77    fdigits(0),
78    ftimeCFD(new TArrayI(24)), 
79    ftimeLED (new TArrayI(24)), 
80    fADC(new TArrayI(24)), 
81    fADC0 (new TArrayI(24)),
82    fSumMult(0),
83    fAmpLED(0),
84    fParam(0)
85 {
86 // ctor which should be used
87
88   AliDebug(1,"processed");
89   fParam = AliT0Parameters::Instance();
90   fParam->Init();
91   Int_t index[25000];
92   Bool_t down=true;
93
94   for (Int_t i=0; i<24; i++){
95     TGraph* gr = fParam ->GetAmpLEDRec(i);
96     Int_t np = gr->GetN();
97     Double_t *x = gr->GetX();
98     Double_t *y = gr->GetY();
99
100     Double_t *x1 = new Double_t[np];
101     Double_t *y1 = new Double_t[np];
102     for (Int_t ii=0; ii<np; ii++) {
103       y1[ii]=y[np-ii-1]; x1[ii]=x[np-ii-1];
104     }
105     
106     TGraph *grInverse = new TGraph(np,y1,x1);
107     fAmpLED.AddAtAndExpand(grInverse,i);
108
109      TGraph* grw = fParam ->GetWalk(i);
110       Int_t npw = grw->GetN();
111      Double_t *yw = grw->GetY();
112     TMath::Sort(npw, yw, index,down);
113     fMaxValue[i]=Int_t(yw[index[0]]);
114   }
115
116     //    delete [] x1;
117     //   delete [] y1;
118 }
119
120 //------------------------------------------------------------------------
121 AliT0Digitizer::~AliT0Digitizer()
122 {
123 // Destructor
124
125   AliDebug(1,"T0");
126
127   delete ftimeCFD;
128   delete fADC;
129   delete ftimeLED;
130   delete  fADC0;
131   //  delete fAmpLED;
132 }
133
134 //------------------------------------------------------------------------
135 Bool_t AliT0Digitizer::Init()
136 {
137 // Initialization
138   AliDebug(1," Init");
139  return kTRUE;
140 }
141  
142 //---------------------------------------------------------------------
143 void AliT0Digitizer::Exec(Option_t* /*option*/)
144 {
145
146   /*
147     Produde digits from hits
148         digits is TObject and includes
149         We are writing array if C & A  TDC
150         C & A  ADC (will need for slow simulation)
151         TOF first particle C & A
152         mean time and time difference (vertex position)
153         
154   */
155
156   //output loader 
157   AliRunLoader *outRL = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
158   AliLoader * pOutStartLoader = outRL->GetLoader("T0Loader");
159
160   AliDebug(1,"start...");
161   //input loader
162   //
163   // From hits to digits
164   //
165   Int_t hit, nhits;
166   Int_t countE[24];
167   Int_t volume, pmt, trCFD, trLED; 
168   Float_t sl, qt;
169   Int_t  bestATDC=0;
170   Int_t  bestCTDC=0;
171   Int_t qtCh=0;
172   Float_t time[24], besttime[24], timeGaus[24] ;
173   //Q->T-> coefficients !!!! should be asked!!!
174   Float_t timeDelayCFD[24];
175   Int_t threshold =50; //photoelectrons
176   Float_t zdetA, zdetC;
177   Int_t sumMultCoeff = 100;
178   Int_t refpoint=0;
179    
180   Int_t ph2Mip = fParam->GetPh2Mip();     
181   Float_t channelWidth = fParam->GetChannelWidth() ;  
182   Float_t delayVertex = fParam->GetTimeDelayTVD();
183    
184   zdetC = TMath::Abs(fParam->GetZPosition("T0/C/PMT1"));
185   zdetA  = TMath::Abs(fParam->GetZPosition("T0/A/PMT15"));
186   
187   AliT0hit  *startHit;
188   TBranch *brHits=0;
189   
190   Int_t nFiles=fManager->GetNinputs();
191   for (Int_t inputFile=0; inputFile<nFiles;  inputFile++) {
192     if (inputFile < nFiles-1) {
193       AliWarning(Form("ignoring input stream %d", inputFile));
194       continue;
195       
196     }
197     
198     Float_t besttimeC=99999.;
199     Float_t besttimeA=99999.;
200     Int_t pmtBestC=99999;
201     Int_t pmtBestA=99999;
202     Int_t timeDiff=99999, meanTime=99999;
203     Int_t sumMult =0;   fSumMult=0;
204     bestATDC = 99999;  bestCTDC = 99999;
205  
206
207     ftimeCFD -> Reset();
208     fADC -> Reset();
209     fADC0 -> Reset();
210     ftimeLED ->Reset();
211     for (Int_t i0=0; i0<24; i0++)
212       {
213         time[i0]=besttime[i0]=timeGaus[i0]=999999; countE[i0]=0;
214       }
215     AliRunLoader * inRL = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
216     AliLoader * pInStartLoader = inRL->GetLoader("T0Loader");
217     if (!inRL->GetAliRun()) inRL->LoadgAlice();
218     fT0  = (AliT0*)inRL ->GetAliRun()->GetDetector("T0");
219
220        //read Hits 
221     pInStartLoader->LoadHits("READ");//probably it is necessary to load them before
222     fHits = fT0->Hits ();
223     TTree *th = pInStartLoader->TreeH();
224     brHits = th->GetBranch("T0");
225     if (brHits) {
226       fT0->SetHitsAddressBranch(brHits);
227     }else{
228       AliWarning("Branch T0 hit not found for this event");
229       // fT0->AddDigit(bestATDC,bestCTDC,meanTime,timeDiff,fSumMult, refpoint,
230       //               ftimeCFD,fADC0,ftimeLED,fADC);
231       continue;      
232     } 
233     Int_t ntracks    = (Int_t) th->GetEntries();
234     
235     if (ntracks<=0) return;
236     // Start loop on tracks in the hits containers
237     for (Int_t track=0; track<ntracks;track++) {
238       brHits->GetEntry(track);
239       nhits = fHits->GetEntriesFast();
240       for (hit=0;hit<nhits;hit++) 
241         {
242           startHit   = (AliT0hit*) fHits->UncheckedAt(hit);
243           if (!startHit) {
244             AliError("The unchecked hit doesn't exist");
245             break;
246           }
247           pmt=startHit->Pmt();
248           Int_t numpmt=pmt-1;
249           Double_t e=startHit->Etot();
250           volume = startHit->Volume();
251           
252           //      if(e>0 && RegisterPhotoE(numpmt,e)) {
253           if(e>0 ) {
254             countE[numpmt]++;
255             besttime[numpmt] = startHit->Time();
256             if(besttime[numpmt]<time[numpmt])
257               {
258                 time[numpmt]=besttime[numpmt];
259               }
260           } //photoelectron accept 
261         } //hits loop
262     } //track loop
263     
264     //spread time A&C by 25ps   && besttime
265     Float_t c = 0.0299792; // cm/ps
266     
267     Float_t koef=(zdetA-zdetC)/c; //correction position difference by cable
268     for (Int_t ipmt=0; ipmt<12; ipmt++){
269       if(countE[ipmt] > threshold) {
270         timeGaus[ipmt]=gRandom->Gaus(time[ipmt],25)+koef;
271         if(timeGaus[ipmt]<besttimeC){
272           besttimeC=timeGaus[ipmt]; //timeC
273           pmtBestC=ipmt;}
274       }
275     }
276     for ( Int_t ipmt=12; ipmt<24; ipmt++){
277       if(countE[ipmt] > threshold) {
278         timeGaus[ipmt]=gRandom->Gaus(time[ipmt],25); 
279         if(timeGaus[ipmt]<besttimeA) {
280           besttimeA=timeGaus[ipmt]; //timeA
281           pmtBestA=ipmt;}
282       } 
283     }
284
285     timeDelayCFD[0] = fParam->GetTimeDelayCFD(0);
286  
287     for (Int_t i=0; i<24; i++)
288       {
289         Float_t  al = countE[i]; 
290         if (al>threshold && timeGaus[i]<50000 ) {
291           //fill ADC
292           // QTC procedure:
293           // phe -> mV 0.3; 1MIP ->500phe -> ln (amp (mV)) = 5;
294           // max 200ns, HIJING  mean 50000phe -> 15000mv -> ln = 15 (s zapasom)
295           // channel 25ps
296           qt= 50.*al/ph2Mip;  // 50mv/Mip amp in mV 
297           //  fill TDC
298           timeDelayCFD[i] = fParam->GetTimeDelayCFD(i);
299           trCFD = Int_t (timeGaus[i]/channelWidth + timeDelayCFD[i]); 
300           TGraph* gr = ((TGraph*)fAmpLED.At(i));
301           sl = gr->Eval(qt);
302
303           trLED = Int_t(( timeGaus[i] + 1000*sl )/channelWidth);
304           qtCh=Int_t (1000.*TMath::Log(qt) / channelWidth);
305           fADC0->AddAt(0,i);
306           fADC->AddAt(qtCh,i);
307           ftimeLED->AddAt(trLED,i); 
308           //      sumMult += Int_t ((al*gain[i]/ph2Mip)*50) ;
309           sumMult += Int_t (qt/sumMultCoeff)  ;
310          
311           // put slewing 
312           TGraph *fu=(TGraph*) fParam ->GetWalk(i) ;
313           Float_t slew=fu->Eval(Float_t(qtCh));
314
315           //      trCFD=trCFD-Int_t(fMaxValue[i]-slew);
316           trCFD = trCFD-Int_t(fMaxValue[i]-slew) + 2000; //for the same channel as cosmic
317           ftimeCFD->AddAt(Int_t (trCFD),i);
318           AliDebug(10,Form("  pmt %i : time in ns %f time in channels %i  LEd %i  ",  i, timeGaus[i],trCFD, trLED ));
319           AliDebug(10,Form(" qt in mV %f qt in ns %f qt in channels %i   ",qt, 
320                            TMath::Log(qt), qtCh));
321
322         }
323       } //pmt loop
324
325     //folding with alignmentz position distribution  
326     if( besttimeC > 10000. && besttimeC <15000)
327       bestCTDC=Int_t ((besttimeC+timeDelayCFD[pmtBestC])
328                          /channelWidth);
329  
330     if( besttimeA > 10000. && besttimeA <15000)
331       bestATDC=Int_t ((besttimeA+timeDelayCFD[pmtBestA])
332                         /channelWidth);
333
334     if (bestATDC < 99999 && bestCTDC < 99999)
335       {
336         timeDiff=Int_t (((besttimeC-besttimeA)+1000*delayVertex)
337                         /channelWidth);
338         meanTime=Int_t (((besttimeC+besttimeA)/2. )/channelWidth);
339       }
340         AliDebug(10,Form(" time A& C %i %i  time diff && mean time in channels %i %i",bestATDC,bestCTDC, timeDiff, meanTime));
341
342     if (sumMult > threshold){
343       fSumMult =  Int_t (1000.* TMath::Log(Double_t(sumMult) / Double_t(sumMultCoeff))
344                          /channelWidth);
345       AliDebug(10,Form("summult mv %i   mult  in chammens %i in ps %i ", 
346                       sumMult, fSumMult, fSumMult*channelWidth));
347     }
348
349     fT0->AddDigit(bestATDC,bestCTDC,meanTime,timeDiff,fSumMult, refpoint,
350                        ftimeCFD,fADC0,ftimeLED,fADC);
351      
352     
353       AliDebug(10,Form(" Digits wrote refpoint %i bestATDC %i bestCTDC %i  meanTime %i  timeDiff %i fSumMult %i ",refpoint ,bestATDC,bestCTDC,meanTime,timeDiff,fSumMult ));
354     pOutStartLoader->UnloadHits();
355   } //input streams loop
356   
357     //load digits    
358     pOutStartLoader->LoadDigits("UPDATE");
359     TTree *treeD  = pOutStartLoader->TreeD();
360     if (treeD == 0x0) {
361       pOutStartLoader->MakeTree("D");
362       treeD = pOutStartLoader->TreeD();
363     }
364     treeD->Reset();
365     fT0  = (AliT0*)outRL ->GetAliRun()->GetDetector("T0");
366     // Make a branch in the tree 
367     fT0->MakeBranch("D");
368      treeD->Fill();
369   
370      pOutStartLoader->WriteDigits("OVERWRITE");
371      
372      fT0->ResetDigits();
373      pOutStartLoader->UnloadDigits();
374      
375 }