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