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