]> git.uio.no Git - u/mrichter/AliRoot.git/blob - START/AliSTARTDigitizer.cxx
Updated comments for Doxygen
[u/mrichter/AliRoot.git] / START / AliSTARTDigitizer.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
18 #include <TTree.h> 
19 #include <TFile.h>
20 #include <TDirectory.h>
21 #include <TRandom.h>
22 #include <TArrayI.h>
23 #include <TError.h>
24 #include <TH1F.h>
25 #include <TGraph.h>
26
27 #include "AliLog.h"
28 #include "AliSTARTDigitizer.h"
29 #include "AliSTART.h"
30 #include "AliSTARThit.h"
31 #include "AliSTARTdigit.h"
32 #include "AliRunDigitizer.h"
33 #include <AliDetector.h>
34 #include "AliRun.h"
35 #include <AliLoader.h>
36 #include <AliRunLoader.h>
37 #include <stdlib.h>
38 #include <Riostream.h>
39 #include <Riostream.h>
40 #include "AliSTARTParameters.h"
41 #include "AliCDBLocal.h"
42 #include "AliCDBStorage.h"
43 #include "AliCDBManager.h"
44 #include "AliCDBEntry.h"
45
46 ClassImp(AliSTARTDigitizer)
47
48 //___________________________________________
49   AliSTARTDigitizer::AliSTARTDigitizer()  :AliDigitizer()
50 {
51 // Default ctor - don't use it
52   ;
53 }
54
55 //___________________________________________
56 AliSTARTDigitizer::AliSTARTDigitizer(AliRunDigitizer* manager) 
57   :AliDigitizer(manager),
58    fSTART(0),
59    fHits(0),
60    fdigits(0),
61    ftimeCFD(0),
62    ftimeLED(0),
63    fADC(0),
64    fADC0(0)
65 {
66 // ctor which should be used
67
68   AliDebug(1,"processed");
69
70   fSTART = 0;
71   fHits = 0;
72   fdigits = 0;
73
74   ftimeCFD = new TArrayI(24); 
75   fADC = new TArrayI(24); 
76   ftimeLED = new TArrayI(24); 
77   fADC0 = new TArrayI(24); 
78   
79
80 }
81
82 //------------------------------------------------------------------------
83 AliSTARTDigitizer::~AliSTARTDigitizer()
84 {
85 // Destructor
86
87   AliDebug(1,"START");
88
89   delete ftimeCFD;
90   delete fADC;
91   delete ftimeLED;
92   delete  fADC0;
93 }
94
95  //------------------------------------------------------------------------
96 Bool_t AliSTARTDigitizer::Init()
97 {
98 // Initialization
99   AliDebug(1," Init");
100  return kTRUE;
101 }
102  
103
104 //---------------------------------------------------------------------
105
106 void AliSTARTDigitizer::Exec(Option_t* /*option*/)
107 {
108
109   /*
110     Produde digits from hits
111         digits is TObject and includes
112         We are writing array if left & right  TDC
113         left & right  ADC (will need for slow simulation)
114         TOF first particle left & right
115         mean time and time difference (vertex position)
116         
117   */
118
119   //output loader 
120   AliRunLoader *outRL = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
121   AliLoader * pOutStartLoader = outRL->GetLoader("STARTLoader");
122
123   AliDebug(1,"start...");
124   //input loader
125   //
126   // From hits to digits
127   //
128  Int_t hit, nhits;
129   Int_t countE[24];
130   Int_t volume, pmt, trCFD, trLED; 
131   Float_t sl, qt;
132   Int_t  bestRightTDC, bestLeftTDC, qtCh;
133   Float_t time[24], besttime[24], timeGaus[24] ;
134     //Q->T-> coefficients !!!! should be asked!!!
135   Float_t gain[24],timeDelayCFD[24], timeDelayLED[24];
136   Int_t threshold =50; //photoelectrons
137   Float_t zdetA, zdetC;
138    Int_t sumMultCoeff = 25;
139   TObjArray slewingLED;
140   TObjArray slewingRec;
141   AliSTARTParameters* param = AliSTARTParameters::Instance();
142   param->Init();
143
144   Int_t ph2Mip = param->GetPh2Mip();     
145   Int_t channelWidth = param->GetChannelWidth() ;  
146   Float_t delayVertex = param->GetTimeDelayTVD();
147   for (Int_t i=0; i<24; i++){
148     timeDelayCFD[i] = param->GetTimeDelayCFD(i);
149     timeDelayLED[i] = param->GetTimeDelayLED(i);
150     gain[i] = param->GetGain(i);
151     TGraph* gr = param ->GetSlew(i);
152     slewingLED.AddAtAndExpand(gr,i);
153
154     TGraph* gr1 = param ->GetSlewRec(i);
155     slewingRec.AddAtAndExpand(gr1,i);
156
157     TGraph* grEff = param ->GetPMTeff(i);
158     fEffPMT.AddAtAndExpand(grEff,i);
159   }
160   zdetC = param->GetZposition(0);
161   zdetA = param->GetZposition(1);
162   
163   AliSTARThit  *startHit;
164   TBranch *brHits=0;
165   
166   Int_t nFiles=fManager->GetNinputs();
167   for (Int_t inputFile=0; inputFile<nFiles;  inputFile++) {
168     if (inputFile < nFiles-1) {
169         AliWarning(Form("ignoring input stream %d", inputFile));
170         continue;
171         
172     }
173     
174     Float_t besttimeright=99999.;
175     Float_t besttimeleft=99999.;
176     Int_t pmtBestRight=9999;
177     Int_t pmtBestLeft=9999;
178     Int_t timeDiff=999, meanTime=0;
179     Int_t sumMult =0;// fSumMult=0;
180     bestRightTDC = 99999;  bestLeftTDC = 99999;
181  
182     ftimeCFD -> Reset();
183     fADC -> Reset();
184     fADC0 -> Reset();
185     ftimeLED ->Reset();
186     for (Int_t i0=0; i0<24; i0++)
187       {
188         time[i0]=besttime[i0]=timeGaus[i0]=99999; countE[i0]=0;
189       }
190     AliRunLoader * inRL = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
191     AliLoader * pInStartLoader = inRL->GetLoader("STARTLoader");
192     if (!inRL->GetAliRun()) inRL->LoadgAlice();
193     AliSTART *fSTART  = (AliSTART*)inRL ->GetAliRun()->GetDetector("START");
194
195        //read Hits 
196     pInStartLoader->LoadHits("READ");//probably it is necessary to load them before
197     TClonesArray *fHits = fSTART->Hits ();
198     TTree *th = pInStartLoader->TreeH();
199     brHits = th->GetBranch("START");
200     if (brHits) {
201       fSTART->SetHitsAddressBranch(brHits);
202     }else{
203       AliError("Branch START hit not found");
204       exit(111);
205     } 
206     Int_t ntracks    = (Int_t) th->GetEntries();
207     
208     if (ntracks<=0) return;
209     // Start loop on tracks in the hits containers
210     for (Int_t track=0; track<ntracks;track++) {
211       brHits->GetEntry(track);
212       nhits = fHits->GetEntriesFast();
213       for (hit=0;hit<nhits;hit++) 
214         {
215           startHit   = (AliSTARThit*) fHits->UncheckedAt(hit);
216           if (!startHit) {
217             AliError("The unchecked hit doesn't exist");
218             break;
219           }
220           pmt=startHit->Pmt();
221           Int_t numpmt=pmt-1;
222                   Double_t e=startHit->Etot();
223           volume = startHit->Volume();
224           
225           if(e>0 && RegisterPhotoE(numpmt,e)) {
226             countE[numpmt]++;
227             besttime[numpmt] = startHit->Time();
228             if(besttime[numpmt]<time[numpmt])
229               {
230                 time[numpmt]=besttime[numpmt];
231               }
232           } //photoelectron accept 
233         } //hits loop
234     } //track loop
235     
236     //spread time right&left by 25ps   && besttime
237     Float_t c = 0.0299792; // cm/ps
238     
239     Float_t koef=(zdetA-zdetC)/c; //correction position difference by cable
240     for (Int_t ipmt=0; ipmt<12; ipmt++){
241       if(countE[ipmt] > threshold) {
242         timeGaus[ipmt]=gRandom->Gaus(time[ipmt],25)+koef;
243         if(timeGaus[ipmt]<besttimeleft){
244           besttimeleft=timeGaus[ipmt]; //timeleft
245           pmtBestLeft=ipmt;}
246      }
247     }
248      for ( Int_t ipmt=12; ipmt<24; ipmt++){
249       if(countE[ipmt] > threshold) {
250         timeGaus[ipmt]=gRandom->Gaus(time[ipmt],25); 
251         if(timeGaus[ipmt]<besttimeright) {
252           besttimeright=timeGaus[ipmt]; //timeright
253           pmtBestRight=ipmt;}
254       } 
255     }
256    //folding with alignmentz position distribution  
257     if( besttimeleft > 10000. && besttimeleft <15000)
258       bestLeftTDC=Int_t ((besttimeleft+1000*timeDelayCFD[pmtBestLeft])
259                          /channelWidth);
260  
261     if( besttimeright > 10000. && besttimeright <15000)
262       bestRightTDC=Int_t ((besttimeright+1000*timeDelayCFD[pmtBestRight])
263                         /channelWidth);
264
265     if (bestRightTDC < 99999 && bestLeftTDC < 99999)
266       {
267         timeDiff=Int_t (((besttimeleft-besttimeright)+1000*delayVertex)
268                         /channelWidth);
269         meanTime=Int_t (((besttimeright+1000*timeDelayCFD[pmtBestLeft]+
270                           besttimeleft+1000*timeDelayCFD[pmtBestLeft])/2.)
271                         /channelWidth);
272         AliDebug(10,Form(" time right& left %i %i  time diff && mean time in channels %i %i",bestRightTDC,bestLeftTDC, timeDiff, meanTime));
273       }
274     for (Int_t i=0; i<24; i++)
275       {
276         Float_t  al = countE[i]; 
277         if (al>threshold && timeGaus[i]<50000 ) {
278           //fill ADC
279           // QTC procedure:
280           // phe -> mV 0.3; 1MIP ->500phe -> ln (amp (mV)) = 5;
281           // max 200ns, HIJING  mean 50000phe -> 15000mv -> ln = 15 (s zapasom)
282           // channel 25ps
283           qt= 50.*al*gain[i]/ph2Mip;  // 50mv/Mip amp in mV 
284           //  fill TDC
285           trCFD = Int_t (timeGaus[i] + 1000.*timeDelayCFD[i])/channelWidth; 
286           trLED= Int_t (timeGaus[i] + 1000.*timeDelayLED[i]); 
287           sl = ((TGraph*)slewingLED.At(i))->Eval(qt);
288           trLED = Int_t(( trLED + 1000*sl )/channelWidth);
289           qtCh=Int_t (1000.*TMath::Log(qt)) / channelWidth;
290           fADC0->AddAt(0,i);
291           fADC->AddAt(qtCh,i);
292           ftimeCFD->AddAt(Int_t (trCFD),i);
293           ftimeLED->AddAt(trLED,i); 
294           //      sumMult += Int_t ((al*gain[i]/ph2Mip)*50) ;
295           sumMult += Int_t (qt/sumMultCoeff)  ;
296          
297         AliDebug(10,Form("  pmt %i : time in ns %f time in channels %i   ",
298                         i, timeGaus[i],trCFD ));
299         AliDebug(10,Form(" qt in mV %f qt in ns %f qt in channels %i   ",qt, 
300                         TMath::Log(qt), qtCh));
301         }
302       } //pmt loop
303
304     if (sumMult > threshold){
305       fSumMult =  Int_t (1000.* TMath::Log(Double_t(sumMult) / Double_t(sumMultCoeff))
306                          /channelWidth);
307       AliDebug(10,Form("summult mv %i   mult  in chammens %i in ps %i ", 
308                       sumMult, fSumMult, fSumMult*channelWidth));
309     }
310     if (  besttimeright<99999 || besttimeleft < 99999) {
311       fSTART->AddDigit(bestRightTDC,bestLeftTDC,meanTime,timeDiff,fSumMult,
312                        ftimeCFD,fADC,ftimeLED,fADC0);
313       
314       AliDebug(10,Form(" Digits wrote bestRightTDC %i bestLeftTDC %i  meanTime %i  timeDiff %i fSumMult %i ", bestRightTDC,bestLeftTDC,meanTime,timeDiff,fSumMult ));
315     }  
316     pOutStartLoader->UnloadHits();
317   } //input streams loop
318   
319     //load digits    
320     pOutStartLoader->LoadDigits("UPDATE");
321     TTree *treeD  = pOutStartLoader->TreeD();
322     if (treeD == 0x0) {
323       pOutStartLoader->MakeTree("D");
324       treeD = pOutStartLoader->TreeD();
325     }
326     treeD->Reset();
327     fSTART  = (AliSTART*)outRL ->GetAliRun()->GetDetector("START");
328     // Make a branch in the tree 
329     fSTART->MakeBranch("D");
330      treeD->Fill();
331   
332      pOutStartLoader->WriteDigits("OVERWRITE");
333      
334      fSTART->ResetDigits();
335      pOutStartLoader->UnloadDigits();
336      
337 }
338
339
340 //------------------------------------------------------------------------
341 Bool_t AliSTARTDigitizer::RegisterPhotoE(Int_t ipmt,Double_t energy)
342 {
343
344   
345   //  Float_t hc=197.326960*1.e6; //mev*nm
346   Double_t hc=1.973*1.e-6; //gev*nm
347   Float_t lambda=hc/energy;
348   Float_t eff = ((TGraph*) fEffPMT.At(ipmt))->Eval(lambda);
349   Double_t  p = gRandom->Rndm();
350
351   if (p > eff)
352     return kFALSE;
353   
354   return kTRUE;
355 }
356
357 //----------------------------------------------------------------------------