]> git.uio.no Git - u/mrichter/AliRoot.git/blob - START/AliSTARTDigitizer.cxx
new digitization and reconstruction corresponded to new data format
[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;
133   Float_t time[24], besttime[24], timeGaus[24] ;
134   //  Float_t channelWidth=25.; //ps
135     //Q->T-> coefficients !!!! should be asked!!!
136   Float_t gain[24],timeDelayCFD[24], timeDelayLED[24];
137   Int_t threshold =50; //photoelectrons
138   Int_t qtThreshold=13, qtCh;
139   Float_t zdetA, zdetC;
140   TObjArray slewingLED;
141   AliSTARTParameters* param = AliSTARTParameters::Instance();
142   Int_t ph2Mip = param->GetPh2Mip();     
143   Int_t channelWidth = param->GetChannelWidth() ;  
144   
145   param->Init();
146   for (Int_t i=0; i<24; i++){
147     timeDelayCFD[i] = param->GetTimeDelayCFD(i);
148     timeDelayLED[i] = param->GetTimeDelayLED(i);
149     gain[i] = param->GetGain(i);
150     TGraph* gr = param ->GetSlew(i);
151     slewingLED.AddAtAndExpand(gr,i);
152     TGraph* grEff = param ->GetPMTeff(i);
153     fEffPMT.AddAtAndExpand(grEff,i);
154   }
155   zdetC = param->GetZposition(0);
156   zdetA = param->GetZposition(1);
157   
158   AliSTARThit  *startHit;
159   TBranch *brHits=0;
160   
161   Int_t nFiles=fManager->GetNinputs();
162   for (Int_t inputFile=0; inputFile<nFiles;  inputFile++) {
163     if (inputFile < nFiles-1) {
164         AliWarning(Form("ignoring input stream %d", inputFile));
165         continue;
166         
167     }
168     
169     Float_t besttimeright=99999.;
170     Float_t besttimeleft=99999.;
171     Int_t pmtBestRight=9999;
172     Int_t pmtBestLeft=9999;
173     Int_t timeDiff=99999, meanTime=0;
174     Int_t sumMult =0;
175     ftimeCFD -> Reset();
176     fADC -> Reset();
177     fADC0 -> Reset();
178     ftimeLED ->Reset();
179     cout<<" before reading event "<<endl;
180     for (Int_t i0=0; i0<24; i0++)
181       {
182         time[i0]=besttime[i0]=timeGaus[i0]=99999; countE[i0]=0;
183       }
184     AliRunLoader * inRL = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(inputFile));
185     AliLoader * pInStartLoader = inRL->GetLoader("STARTLoader");
186     if (!inRL->GetAliRun()) inRL->LoadgAlice();
187     AliSTART *fSTART  = (AliSTART*)inRL ->GetAliRun()->GetDetector("START");
188
189        //read Hits 
190     pInStartLoader->LoadHits("READ");//probably it is necessary to load them before
191     TClonesArray *fHits = fSTART->Hits ();
192     TTree *th = pInStartLoader->TreeH();
193     brHits = th->GetBranch("START");
194     if (brHits) {
195       fSTART->SetHitsAddressBranch(brHits);
196     }else{
197       AliError("Branch START hit not found");
198       exit(111);
199     } 
200     Int_t ntracks    = (Int_t) th->GetEntries();
201     
202     if (ntracks<=0) return;
203     // Start loop on tracks in the hits containers
204     for (Int_t track=0; track<ntracks;track++) {
205       brHits->GetEntry(track);
206       nhits = fHits->GetEntriesFast();
207       for (hit=0;hit<nhits;hit++) 
208         {
209           startHit   = (AliSTARThit*) fHits->UncheckedAt(hit);
210           if (!startHit) {
211             AliError("The unchecked hit doesn't exist");
212             break;
213           }
214           pmt=startHit->Pmt();
215           Int_t numpmt=pmt-1;
216                   Double_t e=startHit->Etot();
217           volume = startHit->Volume();
218           
219           if(e>0 && RegisterPhotoE(numpmt,e)) {
220             countE[numpmt]++;
221             besttime[numpmt] = startHit->Time();
222             if(besttime[numpmt]<time[numpmt])
223               {
224                 time[numpmt]=besttime[numpmt];
225               }
226           } //photoelectron accept 
227         } //hits loop
228     } //track loop
229     cout<<" hits was read "<<endl;
230     
231     //spread time right&left by 25ps   && besttime
232      Float_t c = 0.0299792; // cm/ps
233
234    Float_t koef=(zdetA-zdetC)/c; //correction position difference by cable
235     for (Int_t ipmt=0; ipmt<12; ipmt++){
236       if(countE[ipmt] > threshold) {
237         timeGaus[ipmt]=gRandom->Gaus(time[ipmt],25)+koef;
238         if(timeGaus[ipmt]<besttimeleft){
239           besttimeleft=timeGaus[ipmt]; //timeleft
240           pmtBestLeft=ipmt;}
241       }
242     }
243     cout<<"left "<<besttimeleft<<" "<<pmtBestLeft<<" koef "<<koef<<endl;
244     for ( Int_t ipmt=12; ipmt<24; ipmt++){
245       if(countE[ipmt] > threshold) {
246         timeGaus[ipmt]=gRandom->Gaus(time[ipmt],25); 
247         if(timeGaus[ipmt]<besttimeright) {
248           besttimeright=timeGaus[ipmt]; //timeright
249           pmtBestRight=ipmt;}
250       } 
251     }
252     cout<<"right"<<besttimeright<<" "<<pmtBestRight<<endl;
253   
254    //folding with alignmentz position distribution  
255
256     //    Float_t koef=(zdetA-zdetC)/c;
257     //  Float_t  besttimeleftR= besttimeleft;
258     // besttimeleft=koef+besttimeleft;
259     cout<<"  besttimeleft "<< besttimeleft<<" delay "<< timeDelayCFD[pmtBestLeft]<<" channelWidth "<<channelWidth<<endl;
260     bestLeftTDC=Int_t ((besttimeleft+1000*timeDelayCFD[pmtBestLeft])
261                        /channelWidth);
262     cout<<" bestLeftTDC "<<bestLeftTDC<<" delay "<<timeDelayCFD[pmtBestLeft]<<endl;
263     bestRightTDC=Int_t ((besttimeright+1000*timeDelayCFD[pmtBestLeft])
264                         /channelWidth);
265     timeDiff=Int_t (((besttimeleft+1000*timeDelayCFD[pmtBestLeft])-
266                     (besttimeright+1000*timeDelayCFD[pmtBestLeft])
267                     )/channelWidth);
268     meanTime=Int_t (((besttimeright+1000*timeDelayCFD[pmtBestLeft]+
269                       besttimeleft+1000*timeDelayCFD[pmtBestLeft])/2.)
270                     /channelWidth);
271     cout<<"bestLeftTDC"<<bestLeftTDC<<" bestRightTDC "<< bestRightTDC<<
272       " timeDiff "<<timeDiff<<" meanTime= "<<meanTime<<endl;   
273     AliDebug(2,Form(" time in ns %f ", besttimeleft));
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           cout<<al<<" "<<qt<<" trCFD "<<trCFD<<" sl "<<sl<<" trLED "<<trLED<<
296             " qtThreshold "<<qtThreshold<<"  qtCh "<< qtCh<<endl;
297         }
298       } //pmt loop
299
300         if (sumMult > threshold){
301           fSumMult =  Int_t (TMath::Log(sumMult));
302         }
303     fSTART->AddDigit(bestRightTDC,bestLeftTDC,meanTime,timeDiff,fSumMult,
304                      ftimeCFD,fADC,ftimeLED,fADC0);
305     pOutStartLoader->UnloadHits();
306   } //input streams loop
307   
308     //load digits    
309     pOutStartLoader->LoadDigits("UPDATE");
310     TTree *treeD  = pOutStartLoader->TreeD();
311     if (treeD == 0x0) {
312       pOutStartLoader->MakeTree("D");
313       treeD = pOutStartLoader->TreeD();
314     }
315     treeD->Reset();
316     fSTART  = (AliSTART*)outRL ->GetAliRun()->GetDetector("START");
317     // Make a branch in the tree 
318     fSTART->MakeBranch("D");
319     //     treeD->Print();
320      treeD->Fill();
321      cout<<"   treeD->Fill(); "<<endl;
322   
323      pOutStartLoader->WriteDigits("OVERWRITE");
324      
325      cout<<"  pOutStartLoader->WriteDigits "<<endl;
326      fSTART->ResetDigits();
327      pOutStartLoader->UnloadDigits();
328      
329 }
330
331
332 //------------------------------------------------------------------------
333 Bool_t AliSTARTDigitizer::RegisterPhotoE(Int_t ipmt,Double_t energy)
334 {
335
336   
337   //  Float_t hc=197.326960*1.e6; //mev*nm
338   Double_t hc=1.973*1.e-6; //gev*nm
339   Float_t lambda=hc/energy;
340   Float_t eff = ((TGraph*) fEffPMT.At(ipmt))->Eval(lambda);
341   Double_t  p = gRandom->Rndm();
342
343   if (p > eff)
344     return kFALSE;
345   
346   return kTRUE;
347 }
348
349 //----------------------------------------------------------------------------