]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0Digitizer.cxx
New version of CDB access (Alla). The Look-up-table CDB file has to corrected so...
[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];
143   Int_t threshold =50; //photoelectrons
144   Float_t zdetA, zdetC;
145   Int_t sumMultCoeff = 100;
146   TH1F *hr;
147
148   AliT0Parameters* param = AliT0Parameters::Instance();
149   param->Init();
150
151   
152   Int_t ph2Mip = param->GetPh2Mip();     
153   Int_t channelWidth = param->GetChannelWidth() ;  
154   Float_t delayVertex = param->GetTimeDelayTVD();
155    
156   zdetC = TMath::Abs(param->GetZPosition("T0/C/PMT1"));
157   zdetA  = TMath::Abs(param->GetZPosition("T0/A/PMT15"));
158   
159   AliT0hit  *startHit;
160   TBranch *brHits=0;
161   
162   Int_t nFiles=fManager->GetNinputs();
163   for (Int_t inputFile=0; inputFile<nFiles;  inputFile++) {
164     if (inputFile < nFiles-1) {
165       AliWarning(Form("ignoring input stream %d", inputFile));
166       continue;
167       
168     }
169     
170     Float_t besttimeC=99999.;
171     Float_t besttimeA=99999.;
172     Int_t pmtBestC=9999;
173     Int_t pmtBestA=9999;
174     Int_t timeDiff=999, meanTime=0;
175     Int_t sumMult =0;   fSumMult=0;
176     bestATDC = 99999;  bestCTDC = 99999;
177  
178
179     ftimeCFD -> Reset();
180     fADC -> Reset();
181     fADC0 -> Reset();
182     ftimeLED ->Reset();
183     for (Int_t i0=0; i0<24; i0++)
184       {
185         time[i0]=besttime[i0]=timeGaus[i0]=99999; countE[i0]=0;
186       }
187     AliRunLoader * inRL = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
188     AliLoader * pInStartLoader = inRL->GetLoader("T0Loader");
189     if (!inRL->GetAliRun()) inRL->LoadgAlice();
190     AliT0 *fT0  = (AliT0*)inRL ->GetAliRun()->GetDetector("T0");
191
192        //read Hits 
193     pInStartLoader->LoadHits("READ");//probably it is necessary to load them before
194     TClonesArray *fHits = fT0->Hits ();
195     TTree *th = pInStartLoader->TreeH();
196     brHits = th->GetBranch("T0");
197     if (brHits) {
198       fT0->SetHitsAddressBranch(brHits);
199     }else{
200       AliError("Branch T0 hit not found");
201       exit(111);
202     } 
203     Int_t ntracks    = (Int_t) th->GetEntries();
204     
205     if (ntracks<=0) return;
206     // Start loop on tracks in the hits containers
207     for (Int_t track=0; track<ntracks;track++) {
208       brHits->GetEntry(track);
209       nhits = fHits->GetEntriesFast();
210       for (hit=0;hit<nhits;hit++) 
211         {
212           startHit   = (AliT0hit*) fHits->UncheckedAt(hit);
213           if (!startHit) {
214             AliError("The unchecked hit doesn't exist");
215             break;
216           }
217           pmt=startHit->Pmt();
218           Int_t numpmt=pmt-1;
219           Double_t e=startHit->Etot();
220           volume = startHit->Volume();
221           
222           //      if(e>0 && RegisterPhotoE(numpmt,e)) {
223           if(e>0 ) {
224             countE[numpmt]++;
225             besttime[numpmt] = startHit->Time();
226             if(besttime[numpmt]<time[numpmt])
227               {
228                 time[numpmt]=besttime[numpmt];
229               }
230           } //photoelectron accept 
231         } //hits loop
232     } //track loop
233     
234     //spread time A&C by 25ps   && besttime
235     Float_t c = 0.0299792; // cm/ps
236     
237     Float_t koef=(zdetA-zdetC)/c; //correction position difference by cable
238     for (Int_t ipmt=0; ipmt<12; ipmt++){
239       if(countE[ipmt] > threshold) {
240         timeGaus[ipmt]=gRandom->Gaus(time[ipmt],25)+koef;
241         if(timeGaus[ipmt]<besttimeC){
242           besttimeC=timeGaus[ipmt]; //timeC
243           pmtBestC=ipmt;}
244      }
245     }
246      for ( Int_t ipmt=12; ipmt<24; ipmt++){
247       if(countE[ipmt] > threshold) {
248         timeGaus[ipmt]=gRandom->Gaus(time[ipmt],25); 
249         if(timeGaus[ipmt]<besttimeA) {
250           besttimeA=timeGaus[ipmt]; //timeA
251           pmtBestA=ipmt;}
252       } 
253     }
254    //folding with alignmentz position distribution  
255     if( besttimeC > 10000. && besttimeC <15000)
256       bestCTDC=Int_t ((besttimeC+timeDelayCFD[pmtBestC])
257                          /channelWidth);
258  
259     if( besttimeA > 10000. && besttimeA <15000)
260       bestATDC=Int_t ((besttimeA+timeDelayCFD[pmtBestA])
261                         /channelWidth);
262
263     if (bestATDC < 99999 && bestCTDC < 99999)
264       {
265         timeDiff=Int_t (((besttimeC-besttimeA)+1000*delayVertex)
266                         /channelWidth);
267         meanTime=Int_t (((besttimeC+timeDelayCFD[pmtBestC]+
268                           besttimeA+timeDelayCFD[pmtBestA])/2.)
269                         /channelWidth);
270       }
271         AliDebug(10,Form(" time A& C %i %i  time diff && mean time in channels %i %i",bestATDC,bestCTDC, timeDiff, meanTime));
272           timeDelayCFD[0] = param->GetTimeDelayCFD(0);
273     for (Int_t i=0; i<24; i++)
274       {
275         Float_t  al = countE[i]; 
276         if (al>threshold && timeGaus[i]<50000 ) {
277           //fill ADC
278           // QTC procedure:
279           // phe -> mV 0.3; 1MIP ->500phe -> ln (amp (mV)) = 5;
280           // max 200ns, HIJING  mean 50000phe -> 15000mv -> ln = 15 (s zapasom)
281           // channel 25ps
282           qt= 50.*al/ph2Mip;  // 50mv/Mip amp in mV 
283           //  fill TDC
284           timeDelayCFD[i] = param->GetTimeDelayCFD(i);
285           trCFD = Int_t (timeGaus[i]/channelWidth + (timeDelayCFD[i]-timeDelayCFD[0])); 
286           //      trLED= Int_t (timeGaus[i] + timeDelayLED[i]); 
287           TGraph* gr = param ->GetAmpLED(i);
288           sl = gr->Eval(qt);
289           trLED = Int_t(( timeGaus[i] + 1000*sl )/channelWidth);
290           qtCh=Int_t (1000.*TMath::Log(qt)) / channelWidth;
291           fADC0->AddAt(0,i);
292           fADC->AddAt(qtCh,i);
293           ftimeLED->AddAt(trLED,i); 
294           //      sumMult += Int_t ((al*gain[i]/ph2Mip)*50) ;
295           sumMult += Int_t (qt/sumMultCoeff)  ;
296          
297           // put slewing 
298           TGraph* fu = param ->GetWalk(i);
299           Float_t slew=fu->Eval(Float_t(qtCh));
300           hr=fu->GetHistogram();
301           Float_t maxValue=hr->GetMaximum(50);
302           trCFD=trCFD-Int_t((maxValue-slew)/channelWidth);
303           ftimeCFD->AddAt(Int_t (trCFD),i);
304           AliDebug(10,Form("  pmt %i : time in ns %f time in channels %i   ",
305                            i, timeGaus[i],trCFD ));
306           AliDebug(10,Form(" qt in mV %f qt in ns %f qt in channels %i   ",qt, 
307                            TMath::Log(qt), qtCh));
308
309         }
310       } //pmt loop
311
312     if (sumMult > threshold){
313       fSumMult =  Int_t (1000.* TMath::Log(Double_t(sumMult) / Double_t(sumMultCoeff))
314                          /channelWidth);
315       AliDebug(10,Form("summult mv %i   mult  in chammens %i in ps %i ", 
316                       sumMult, fSumMult, fSumMult*channelWidth));
317     }
318
319       fT0->AddDigit(bestATDC,bestCTDC,meanTime,timeDiff,fSumMult,
320                        ftimeCFD,fADC0,ftimeLED,fADC);
321      
322       AliDebug(10,Form(" Digits wrote bestATDC %i bestCTDC %i  meanTime %i  timeDiff %i fSumMult %i ", bestATDC,bestCTDC,meanTime,timeDiff,fSumMult ));
323     pOutStartLoader->UnloadHits();
324   } //input streams loop
325   
326     //load digits    
327     pOutStartLoader->LoadDigits("UPDATE");
328     TTree *treeD  = pOutStartLoader->TreeD();
329     if (treeD == 0x0) {
330       pOutStartLoader->MakeTree("D");
331       treeD = pOutStartLoader->TreeD();
332     }
333     treeD->Reset();
334     fT0  = (AliT0*)outRL ->GetAliRun()->GetDetector("T0");
335     // Make a branch in the tree 
336     fT0->MakeBranch("D");
337      treeD->Fill();
338   
339      pOutStartLoader->WriteDigits("OVERWRITE");
340      
341      fT0->ResetDigits();
342      pOutStartLoader->UnloadDigits();
343      
344 }