]> git.uio.no Git - u/mrichter/AliRoot.git/blob - T0/AliT0Digitizer.cxx
violations fixed
[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
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();
97
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];
102     }
103     
104     TGraph *grInverse = new TGraph(np,y1,x1);
105     delete [] x1;
106     delete [] y1;
107     fAmpLED.AddAtAndExpand(grInverse,i);
108   }
109 }
110 //------------------------------------------------------------------------
111 AliT0Digitizer::~AliT0Digitizer()
112 {
113 // Destructor
114
115   AliDebug(1,"T0");
116
117   delete ftimeCFD;
118   delete fADC;
119   delete ftimeLED;
120   delete  fADC0;
121   //  delete fAmpLED;
122 }
123
124 //------------------------------------------------------------------------
125 Bool_t AliT0Digitizer::Init()
126 {
127 // Initialization
128   AliDebug(1," Init");
129  return kTRUE;
130 }
131  
132 //---------------------------------------------------------------------
133 void AliT0Digitizer::Exec(Option_t* /*option*/)
134 {
135
136   /*
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)
143         
144   */
145
146   //output loader 
147   AliRunLoader *outRL = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
148   AliLoader * pOutStartLoader = outRL->GetLoader("T0Loader");
149
150   AliDebug(1,"start...");
151   //input loader
152   //
153   // From hits to digits
154   //
155   Int_t hit, nhits;
156   Int_t countE[24];
157   Int_t volume, pmt, trCFD, trLED; 
158   Float_t sl, qt;
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;
166   Int_t refpoint=0;
167   TH1F *hr;
168
169
170   
171   Int_t ph2Mip = fParam->GetPh2Mip();     
172   Float_t channelWidth = fParam->GetChannelWidth() ;  
173   Float_t delayVertex = fParam->GetTimeDelayTVD();
174    
175   zdetC = TMath::Abs(fParam->GetZPosition("T0/C/PMT1"));
176   zdetA  = TMath::Abs(fParam->GetZPosition("T0/A/PMT15"));
177   
178   AliT0hit  *startHit;
179   TBranch *brHits=0;
180   
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));
185       continue;
186       
187     }
188     
189     Float_t besttimeC=99999.;
190     Float_t besttimeA=99999.;
191     Int_t pmtBestC=9999;
192     Int_t pmtBestA=9999;
193     Int_t timeDiff=999, meanTime=0;
194     Int_t sumMult =0;   fSumMult=0;
195     bestATDC = 99999;  bestCTDC = 99999;
196  
197
198     ftimeCFD -> Reset();
199     fADC -> Reset();
200     fADC0 -> Reset();
201     ftimeLED ->Reset();
202     for (Int_t i0=0; i0<24; i0++)
203       {
204         time[i0]=besttime[i0]=timeGaus[i0]=99999; countE[i0]=0;
205       }
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");
210
211        //read Hits 
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");
216     if (brHits) {
217       fT0->SetHitsAddressBranch(brHits);
218     }else{
219       AliError("Branch T0 hit not found");
220       exit(111);
221     } 
222     Int_t ntracks    = (Int_t) th->GetEntries();
223     
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++) 
230         {
231           startHit   = (AliT0hit*) fHits->UncheckedAt(hit);
232           if (!startHit) {
233             AliError("The unchecked hit doesn't exist");
234             break;
235           }
236           pmt=startHit->Pmt();
237           Int_t numpmt=pmt-1;
238           Double_t e=startHit->Etot();
239           volume = startHit->Volume();
240           
241           //      if(e>0 && RegisterPhotoE(numpmt,e)) {
242           if(e>0 ) {
243             countE[numpmt]++;
244             besttime[numpmt] = startHit->Time();
245             if(besttime[numpmt]<time[numpmt])
246               {
247                 time[numpmt]=besttime[numpmt];
248               }
249           } //photoelectron accept 
250         } //hits loop
251     } //track loop
252     
253     //spread time A&C by 25ps   && besttime
254     Float_t c = 0.0299792; // cm/ps
255     
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
262           pmtBestC=ipmt;}
263       }
264     }
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
270           pmtBestA=ipmt;}
271       } 
272     }
273
274     timeDelayCFD[0] = fParam->GetTimeDelayCFD(0);
275     Int_t meanTimeDelay=200;
276
277     for (Int_t i=0; i<24; i++)
278       {
279         Float_t  al = countE[i]; 
280         if (al>threshold && timeGaus[i]<50000 ) {
281           //fill ADC
282           // QTC procedure:
283           // phe -> mV 0.3; 1MIP ->500phe -> ln (amp (mV)) = 5;
284           // max 200ns, HIJING  mean 50000phe -> 15000mv -> ln = 15 (s zapasom)
285           // channel 25ps
286           qt= 50.*al/ph2Mip;  // 50mv/Mip amp in mV 
287           //  fill TDC
288           timeDelayCFD[i] = fParam->GetTimeDelayCFD(i);
289           trCFD = Int_t (timeGaus[i]/channelWidth + (timeDelayCFD[i]-timeDelayCFD[0])); 
290           TGraph* gr = ((TGraph*)fAmpLED.At(i));
291           sl = gr->Eval(qt);
292
293           trLED = Int_t(( timeGaus[i] + 1000*sl )/channelWidth);
294           qtCh=Int_t (1000.*TMath::Log(qt) / channelWidth);
295           fADC0->AddAt(0,i);
296           fADC->AddAt(qtCh,i);
297           ftimeLED->AddAt(trLED,i); 
298           //      sumMult += Int_t ((al*gain[i]/ph2Mip)*50) ;
299           sumMult += Int_t (qt/sumMultCoeff)  ;
300          
301           // put slewing 
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));
312
313         }
314       } //pmt loop
315
316     //folding with alignmentz position distribution  
317     if( besttimeC > 10000. && besttimeC <15000)
318       bestCTDC=Int_t ((besttimeC+timeDelayCFD[pmtBestC])
319                          /channelWidth);
320  
321     if( besttimeA > 10000. && besttimeA <15000)
322       bestATDC=Int_t ((besttimeA+timeDelayCFD[pmtBestA])
323                         /channelWidth);
324
325     if (bestATDC < 99999 && bestCTDC < 99999)
326       {
327         timeDiff=Int_t (((besttimeC-besttimeA)+1000*delayVertex)
328                         /channelWidth);
329         meanTime=Int_t (((besttimeC+besttimeA)/2. + meanTimeDelay)/channelWidth);
330       }
331         AliDebug(10,Form(" time A& C %i %i  time diff && mean time in channels %i %i",bestATDC,bestCTDC, timeDiff, meanTime));
332
333     if (sumMult > threshold){
334       fSumMult =  Int_t (1000.* TMath::Log(Double_t(sumMult) / Double_t(sumMultCoeff))
335                          /channelWidth);
336       AliDebug(10,Form("summult mv %i   mult  in chammens %i in ps %i ", 
337                       sumMult, fSumMult, fSumMult*channelWidth));
338     }
339
340     fT0->AddDigit(bestATDC,bestCTDC,meanTime,timeDiff,fSumMult, refpoint,
341                        ftimeCFD,fADC0,ftimeLED,fADC);
342      
343     
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
347   
348     //load digits    
349     pOutStartLoader->LoadDigits("UPDATE");
350     TTree *treeD  = pOutStartLoader->TreeD();
351     if (treeD == 0x0) {
352       pOutStartLoader->MakeTree("D");
353       treeD = pOutStartLoader->TreeD();
354     }
355     treeD->Reset();
356     fT0  = (AliT0*)outRL ->GetAliRun()->GetDetector("T0");
357     // Make a branch in the tree 
358     fT0->MakeBranch("D");
359      treeD->Fill();
360   
361      pOutStartLoader->WriteDigits("OVERWRITE");
362      
363      fT0->ResetDigits();
364      pOutStartLoader->UnloadDigits();
365      
366 }