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