1 /**************************************************************************
2 * Copyright(c) 1998-2000, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
19 Class for creating of the sumable digits and digits from MC data
21 The input : ideal signals (Hits->Diffusion->Attachment -Ideal signal)
22 The output: raw digits
25 1. Pad by pad gain map
27 3. The dead channels identified - zerro noise for corresponding pads
28 In this case the outpu equal zerro
37 #include <TObjArray.h>
39 #include <TDirectory.h>
40 #include <Riostream.h>
41 #include <TParameter.h>
43 #include "AliTPCDigitizer.h"
46 #include "AliTPCParam.h"
47 #include "AliTPCParamSR.h"
50 #include "AliRunDigitizer.h"
51 #include "AliSimDigits.h"
54 #include "AliTPCcalibDB.h"
55 #include "AliTPCCalPad.h"
56 #include "AliTPCCalROC.h"
58 ClassImp(AliTPCDigitizer)
60 //___________________________________________
61 AliTPCDigitizer::AliTPCDigitizer() :AliDigitizer(),fDebug(0)
64 // Default ctor - don't use it
69 //___________________________________________
70 AliTPCDigitizer::AliTPCDigitizer(AliRunDigitizer* manager)
71 :AliDigitizer(manager),fDebug(0)
74 // ctor which should be used
76 AliDebug(2,"(AliRunDigitizer* manager) was processed");
79 //------------------------------------------------------------------------
80 AliTPCDigitizer::~AliTPCDigitizer()
87 //------------------------------------------------------------------------
88 Bool_t AliTPCDigitizer::Init()
96 //------------------------------------------------------------------------
97 void AliTPCDigitizer::Exec(Option_t* option)
101 //------------------------------------------------------------------------
102 void AliTPCDigitizer::ExecFast(Option_t* option)
105 // merge input tree's with summable digits
106 //output stored in TreeTPCD
109 TString optionString = option;
110 if (!strcmp(optionString.Data(),"deb")) {
111 cout<<"AliTPCDigitizer::Exec: called with option deb "<<endl;
114 //get detector and geometry
117 AliRunLoader *rl, *orl;
118 AliLoader *gime, *ogime;
122 Warning("ExecFast","gAlice is NULL. Loading from input 0");
123 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
126 Error("ExecFast","Can not find Run Loader for input 0. Can not proceed.");
132 AliTPC *pTPC = (AliTPC *) gAlice->GetModule("TPC");
133 AliTPCParam * param = pTPC->GetParam();
135 sprintf(s,param->GetTitle());
136 sprintf(ss,"75x40_100x60");
138 printf("2 pad-length geom hits with 3 pad-lenght geom digits...\n");
140 param=new AliTPCParamSR();
143 sprintf(ss,"75x40_100x60_150x60");
144 if(strcmp(s,ss)!=0) {
145 printf("No TPC parameters found...\n");
150 pTPC->GenerNoise(500000); //create table with noise
152 Int_t nInputs = fManager->GetNinputs();
153 Int_t * masks = new Int_t[nInputs];
154 for (Int_t i=0; i<nInputs;i++)
155 masks[i]= fManager->GetMask(i);
156 Short_t **pdig= new Short_t*[nInputs]; //pointers to the expanded digits array
157 Int_t **ptr= new Int_t*[nInputs]; //pointers to the expanded tracks array
158 Bool_t *active= new Bool_t[nInputs]; //flag for active input segments
161 //create digits array for given sectors
164 //create branch's in TPC treeD
165 orl = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
166 ogime = orl->GetLoader("TPCLoader");
167 TTree * tree = ogime->TreeD();
168 AliSimDigits * digrow = new AliSimDigits;
172 ogime->MakeTree("D");
173 tree = ogime->TreeD();
175 tree->Branch("Segment","AliSimDigits",&digrow);
177 AliSimDigits ** digarr = new AliSimDigits*[nInputs];
178 for (Int_t i1=0;i1<nInputs; i1++)
182 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(i1));
183 gime = rl->GetLoader("TPCLoader");
184 gime->LoadSDigits("read");
185 TTree * treear = gime->TreeS();
189 cerr<<"AliTPCDigitizer: Input tree with SDigits not found in"
190 <<" input "<< i1<<endl;
194 sprintf(phname,"lhcphase%d",i1);
195 TParameter<float> *ph = (TParameter<float>*)treear->GetUserInfo()
196 ->FindObject("lhcphase0");
198 cerr<<"AliTPCDigitizer: LHC phase not found in"
199 <<" input "<< i1<<endl;
202 tree->GetUserInfo()->Add(new TParameter<float>(phname,ph->GetVal()));
204 if (treear->GetIndex()==0)
205 treear->BuildIndex("fSegmentID","fSegmentID");
206 treear->GetBranch("Segment")->SetAddress(&digarr[i1]);
214 param->SetZeroSup(2);
216 Int_t zerosup = param->GetZeroSup();
217 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetDedxGainFactor();
218 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();
220 //Loop over segments of the TPC
222 for (Int_t segmentID=0; segmentID<param->GetNRowsTotal(); segmentID++)
225 if (!param->AdjustSectorRow(segmentID,sec,row))
227 cerr<<"AliTPC warning: invalid segment ID ! "<<segmentID<<endl;
230 AliTPCCalROC * gainROC = gainTPC->GetCalROC(sec); // pad gains per given sector
231 AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(sec); // noise per given sector
232 digrow->SetID(segmentID);
237 Bool_t digitize = kFALSE;
238 for (Int_t i=0;i<nInputs; i++)
241 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(i));
242 gime = rl->GetLoader("TPCLoader");
244 if (gime->TreeS()->GetEntryWithIndex(segmentID,segmentID) >= 0) {
245 digarr[i]->ExpandBuffer();
246 digarr[i]->ExpandTrackBuffer();
247 nrows = digarr[i]->GetNRows();
248 ncols = digarr[i]->GetNCols();
250 if (!fRegionOfInterest || (i == 0)) digitize = kTRUE;
254 if (fRegionOfInterest && !digitize) break;
256 if (!digitize) continue;
258 digrow->Allocate(nrows,ncols);
259 digrow->AllocateTrack(3);
262 Int_t label[1000]; //stack for 300 events
265 Int_t nElems = nrows*ncols;
267 for (Int_t i=0;i<nInputs; i++)
269 pdig[i] = digarr[i]->GetDigits();
270 ptr[i] = digarr[i]->GetTracks();
273 Short_t *pdig1= digrow->GetDigits();
274 Int_t *ptr1= digrow->GetTracks() ;
278 for (Int_t elem=0;elem<nElems; elem++)
284 for (Int_t i=0;i<nInputs; i++) if (active[i])
286 // q += digarr[i]->GetDigitFast(rows,col);
289 for (Int_t tr=0;tr<3;tr++)
291 // Int_t lab = digarr[i]->GetTrackIDFast(rows,col,tr);
292 Int_t lab = ptr[i][tr*nElems];
293 if ( (lab > 1) && *(pdig[i])>zerosup)
295 label[labptr]=lab+masks[i];
302 q/=16.; //conversion factor
303 Float_t gain = gainROC->GetValue(row,elem/nrows); // get gain for given - pad-row pad
305 //printf("problem\n");
308 Float_t noisePad = noiseROC->GetValue(row,elem/nrows);
309 // Float_t noise = gRandom->Gaus(0,param->GetNoise()*param->GetNoiseNormFac());
310 Float_t noise = pTPC->GetNoise();
315 if(q >= param->GetADCSat()) q = (Short_t)(param->GetADCSat() - 1);
316 //digrow->SetDigitFast((Short_t)q,rows,col);
318 for (Int_t tr=0;tr<3;tr++)
321 ptr1[tr*nElems] = label[tr];
330 digrow->GlitchFilter();
332 digrow->CompresBuffer(1,zerosup);
333 digrow->CompresTrackBuffer(1);
335 if (fDebug>0) cerr<<sec<<"\t"<<row<<"\n";
336 } //for (Int_t n=0; n<param->GetNRowsTotal(); n++)
339 orl = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
340 ogime = orl->GetLoader("TPCLoader");
341 ogime->WriteDigits("OVERWRITE");
343 //fManager->GetTreeDTPC()->Write(0,TObject::kOverwrite);
346 for (Int_t i1=0;i1<nInputs; i1++) delete digarr[i1];
356 //------------------------------------------------------------------------
357 void AliTPCDigitizer::ExecSave(Option_t* option)
360 // merge input tree's with summable digits
361 //output stored in TreeTPCD
363 TString optionString = option;
364 if (!strcmp(optionString.Data(),"deb")) {
365 cout<<"AliTPCDigitizer::Exec: called with option deb "<<endl;
368 //get detector and geometry
369 AliRunLoader *rl, *orl;
370 AliLoader *gime, *ogime;
373 orl = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
374 ogime = orl->GetLoader("TPCLoader");
376 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
377 gime = rl->GetLoader("TPCLoader");
380 AliRun* alirun = rl->GetAliRun();
382 AliTPC *pTPC = (AliTPC *) alirun->GetModule("TPC");
383 AliTPCParam * param = pTPC->GetParam();
384 pTPC->GenerNoise(500000); //create teble with noise
385 printf("noise %f \n", param->GetNoise()*param->GetNoiseNormFac());
387 Int_t nInputs = fManager->GetNinputs();
388 Int_t * masks = new Int_t[nInputs];
389 for (Int_t i=0; i<nInputs;i++)
390 masks[i]= fManager->GetMask(i);
392 AliSimDigits ** digarr = new AliSimDigits*[nInputs];
393 for (Int_t i1=0;i1<nInputs; i1++)
397 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(i1));
398 gime = rl->GetLoader("TPCLoader");
400 TTree * treear = gime->TreeS();
401 TBranch * br = treear->GetBranch("fSegmentID");
402 if (br) br->GetFile()->cd();
404 cerr<<" TPC - not existing input = \n"<<i1<<" ";
406 treear->GetBranch("Segment")->SetAddress(&digarr[i1]);
409 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
410 gime = rl->GetLoader("TPCLoader");
411 Stat_t nentries = gime->TreeS()->GetEntries();
414 //create branch's in TPC treeD
415 AliSimDigits * digrow = new AliSimDigits;
416 TTree * tree = ogime->TreeD();
418 tree->Branch("Segment","AliSimDigits",&digrow);
419 param->SetZeroSup(2);
421 Int_t zerosup = param->GetZeroSup();
422 //Loop over segments of the TPC
424 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetDedxGainFactor();
425 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();
426 for (Int_t n=0; n<nentries; n++) {
427 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
428 gime = rl->GetLoader("TPCLoader");
429 gime->TreeS()->GetEvent(n);
431 digarr[0]->ExpandBuffer();
432 digarr[0]->ExpandTrackBuffer();
435 for (Int_t i=1;i<nInputs; i++){
436 // fManager->GetInputTreeTPCS(i)->GetEntryWithIndex(digarr[0]->GetID(),digarr[0]->GetID());
437 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(i));
438 gime = rl->GetLoader("TPCLoader");
439 gime->TreeS()->GetEntryWithIndex(digarr[0]->GetID(),digarr[0]->GetID());
440 digarr[i]->ExpandBuffer();
441 digarr[i]->ExpandTrackBuffer();
442 if ((digarr[0]->GetID()-digarr[i]->GetID())>0)
448 if (!param->AdjustSectorRow(digarr[0]->GetID(),sec,row)) {
449 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr[0]->GetID()<<endl;
453 AliTPCCalROC * gainROC = gainTPC->GetCalROC(sec); // pad gains per given sector
454 AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(sec); // noise per given sector
455 digrow->SetID(digarr[0]->GetID());
457 Int_t nrows = digarr[0]->GetNRows();
458 Int_t ncols = digarr[0]->GetNCols();
459 digrow->Allocate(nrows,ncols);
460 digrow->AllocateTrack(3);
463 Int_t label[1000]; //stack for 300 events
468 for (Int_t rows=0;rows<nrows; rows++){
469 for (Int_t col=0;col<ncols; col++){
474 for (Int_t i=0;i<nInputs; i++){
475 q += digarr[i]->GetDigitFast(rows,col);
478 for (Int_t tr=0;tr<3;tr++) {
479 Int_t lab = digarr[i]->GetTrackIDFast(rows,col,tr);
480 //Int_t lab = ptr[i][tr*nElems];
482 label[labptr]=lab+masks[i];
490 q/=16.; //conversion factor
491 // Float_t noise = gRandom->Gaus(0,param->GetNoise()*param->GetNoiseNormFac());
492 Float_t gain = gainROC->GetValue(row,col);
494 Float_t noisePad = noiseROC->GetValue(row, col);
496 Float_t noise = pTPC->GetNoise();
502 if(q >= param->GetADCSat()) q = (Short_t)(param->GetADCSat() - 1);
503 digrow->SetDigitFast((Short_t)q,rows,col);
504 // *pdig1 =Short_t(q);
505 for (Int_t tr=0;tr<3;tr++){
507 ((AliSimDigits*)digrow)->SetTrackIDFast(label[tr],rows,col,tr);
508 //ptr1[tr*nElems] = label[tr];
510 // ((AliSimDigits*)digrow)->SetTrackIDFast(-1,rows,col,tr);
511 // ptr1[tr*nElems] = 1;
519 digrow->CompresBuffer(1,zerosup);
520 digrow->CompresTrackBuffer(1);
522 if (fDebug>0) cerr<<sec<<"\t"<<row<<"\n";
524 // printf("end TPC merging - end -Tree %s\t%p\n",fManager->GetInputTreeH(0)->GetName(),fManager->GetInputTreeH(0)->GetListOfBranches()->At(3));
525 //fManager->GetTreeDTPC()->Write(0,TObject::kOverwrite);
526 ogime->WriteDigits("OVERWRITE");
528 for (Int_t i=1;i<nInputs; i++)
530 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(i));
531 gime = rl->GetLoader("TPCLoader");
532 gime->UnloadSDigits();
534 ogime->UnloadDigits();
537 for (Int_t i1=0;i1<nInputs; i1++) delete digarr[i1];