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