2 /**************************************************************************
3 * Copyright(c) 1998-2000, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
20 #include <TDirectory.h>
28 #include "AliSTARTDigitizer.h"
30 #include "AliSTARThit.h"
31 #include "AliSTARTdigit.h"
32 #include "AliRunDigitizer.h"
33 #include <AliDetector.h>
35 #include <AliLoader.h>
36 #include <AliRunLoader.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"
46 ClassImp(AliSTARTDigitizer)
48 //___________________________________________
49 AliSTARTDigitizer::AliSTARTDigitizer() :AliDigitizer()
51 // Default ctor - don't use it
55 //___________________________________________
56 AliSTARTDigitizer::AliSTARTDigitizer(AliRunDigitizer* manager)
57 :AliDigitizer(manager),
66 // ctor which should be used
68 AliDebug(1,"processed");
74 ftimeCFD = new TArrayI(24);
75 fADC = new TArrayI(24);
76 ftimeLED = new TArrayI(24);
77 fADC0 = new TArrayI(24);
82 //------------------------------------------------------------------------
83 AliSTARTDigitizer::~AliSTARTDigitizer()
95 //------------------------------------------------------------------------
96 Bool_t AliSTARTDigitizer::Init()
104 //---------------------------------------------------------------------
106 void AliSTARTDigitizer::Exec(Option_t* /*option*/)
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)
120 AliRunLoader *outRL = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
121 AliLoader * pOutStartLoader = outRL->GetLoader("STARTLoader");
123 AliDebug(1,"start...");
126 // From hits to digits
130 Int_t volume, pmt, trCFD, trLED;
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();
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);
154 TGraph* gr1 = param ->GetSlewRec(i);
155 slewingRec.AddAtAndExpand(gr1,i);
157 TGraph* grEff = param ->GetPMTeff(i);
158 fEffPMT.AddAtAndExpand(grEff,i);
160 zdetC = param->GetZposition(0);
161 zdetA = param->GetZposition(1);
163 AliSTARThit *startHit;
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));
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;
186 for (Int_t i0=0; i0<24; i0++)
188 time[i0]=besttime[i0]=timeGaus[i0]=99999; countE[i0]=0;
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");
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");
201 fSTART->SetHitsAddressBranch(brHits);
203 AliError("Branch START hit not found");
206 Int_t ntracks = (Int_t) th->GetEntries();
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++)
215 startHit = (AliSTARThit*) fHits->UncheckedAt(hit);
217 AliError("The unchecked hit doesn't exist");
222 Double_t e=startHit->Etot();
223 volume = startHit->Volume();
225 if(e>0 && RegisterPhotoE(numpmt,e)) {
227 besttime[numpmt] = startHit->Time();
228 if(besttime[numpmt]<time[numpmt])
230 time[numpmt]=besttime[numpmt];
232 } //photoelectron accept
236 //spread time right&left by 25ps && besttime
237 Float_t c = 0.0299792; // cm/ps
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
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
256 //folding with alignmentz position distribution
257 if( besttimeleft > 10000. && besttimeleft <15000)
258 bestLeftTDC=Int_t ((besttimeleft+1000*timeDelayCFD[pmtBestLeft])
261 if( besttimeright > 10000. && besttimeright <15000)
262 bestRightTDC=Int_t ((besttimeright+1000*timeDelayCFD[pmtBestRight])
265 if (bestRightTDC < 99999 && bestLeftTDC < 99999)
267 timeDiff=Int_t (((besttimeleft-besttimeright)+1000*delayVertex)
269 meanTime=Int_t (((besttimeright+1000*timeDelayCFD[pmtBestLeft]+
270 besttimeleft+1000*timeDelayCFD[pmtBestLeft])/2.)
272 AliDebug(10,Form(" time right& left %i %i time diff && mean time in channels %i %i",bestRightTDC,bestLeftTDC, timeDiff, meanTime));
274 for (Int_t i=0; i<24; i++)
276 Float_t al = countE[i];
277 if (al>threshold && timeGaus[i]<50000 ) {
280 // phe -> mV 0.3; 1MIP ->500phe -> ln (amp (mV)) = 5;
281 // max 200ns, HIJING mean 50000phe -> 15000mv -> ln = 15 (s zapasom)
283 qt= 50.*al*gain[i]/ph2Mip; // 50mv/Mip amp in mV
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;
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) ;
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));
304 if (sumMult > threshold){
305 fSumMult = Int_t (1000.* TMath::Log(Double_t(sumMult) / Double_t(sumMultCoeff))
307 AliDebug(10,Form("summult mv %i mult in chammens %i in ps %i ",
308 sumMult, fSumMult, fSumMult*channelWidth));
310 if ( besttimeright<99999 || besttimeleft < 99999) {
311 fSTART->AddDigit(bestRightTDC,bestLeftTDC,meanTime,timeDiff,fSumMult,
312 ftimeCFD,fADC,ftimeLED,fADC0);
314 AliDebug(10,Form(" Digits wrote bestRightTDC %i bestLeftTDC %i meanTime %i timeDiff %i fSumMult %i ", bestRightTDC,bestLeftTDC,meanTime,timeDiff,fSumMult ));
316 pOutStartLoader->UnloadHits();
317 } //input streams loop
320 pOutStartLoader->LoadDigits("UPDATE");
321 TTree *treeD = pOutStartLoader->TreeD();
323 pOutStartLoader->MakeTree("D");
324 treeD = pOutStartLoader->TreeD();
327 fSTART = (AliSTART*)outRL ->GetAliRun()->GetDetector("START");
328 // Make a branch in the tree
329 fSTART->MakeBranch("D");
332 pOutStartLoader->WriteDigits("OVERWRITE");
334 fSTART->ResetDigits();
335 pOutStartLoader->UnloadDigits();
340 //------------------------------------------------------------------------
341 Bool_t AliSTARTDigitizer::RegisterPhotoE(Int_t ipmt,Double_t energy)
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();
357 //----------------------------------------------------------------------------