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