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