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>
42 #include "AliTPCDigitizer.h"
45 #include "AliTPCParam.h"
46 #include "AliTPCParamSR.h"
49 #include "AliRunDigitizer.h"
50 #include "AliSimDigits.h"
53 #include "AliTPCcalibDB.h"
54 #include "AliTPCCalPad.h"
55 #include "AliTPCCalROC.h"
57 ClassImp(AliTPCDigitizer)
59 //___________________________________________
60 AliTPCDigitizer::AliTPCDigitizer() :AliDigitizer(),fDebug(0)
63 // Default ctor - don't use it
68 //___________________________________________
69 AliTPCDigitizer::AliTPCDigitizer(AliRunDigitizer* manager)
70 :AliDigitizer(manager),fDebug(0)
73 // ctor which should be used
75 AliDebug(2,"(AliRunDigitizer* manager) was processed");
78 //------------------------------------------------------------------------
79 AliTPCDigitizer::~AliTPCDigitizer()
86 //------------------------------------------------------------------------
87 Bool_t AliTPCDigitizer::Init()
95 //------------------------------------------------------------------------
96 void AliTPCDigitizer::Exec(Option_t* option)
100 //------------------------------------------------------------------------
101 void AliTPCDigitizer::ExecFast(Option_t* option)
104 // merge input tree's with summable digits
105 //output stored in TreeTPCD
108 TString optionString = option;
109 if (!strcmp(optionString.Data(),"deb")) {
110 cout<<"AliTPCDigitizer::Exec: called with option deb "<<endl;
113 //get detector and geometry
116 AliRunLoader *rl, *orl;
117 AliLoader *gime, *ogime;
121 Warning("ExecFast","gAlice is NULL. Loading from input 0");
122 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
125 Error("ExecFast","Can not find Run Loader for input 0. Can not proceed.");
131 AliTPC *pTPC = (AliTPC *) gAlice->GetModule("TPC");
132 AliTPCParam * param = pTPC->GetParam();
134 sprintf(s,param->GetTitle());
135 sprintf(ss,"75x40_100x60");
137 printf("2 pad-length geom hits with 3 pad-lenght geom digits...\n");
139 param=new AliTPCParamSR();
142 sprintf(ss,"75x40_100x60_150x60");
143 if(strcmp(s,ss)!=0) {
144 printf("No TPC parameters found...\n");
149 pTPC->GenerNoise(500000); //create teble with noise
151 Int_t nInputs = fManager->GetNinputs();
152 Int_t * masks = new Int_t[nInputs];
153 for (Int_t i=0; i<nInputs;i++)
154 masks[i]= fManager->GetMask(i);
155 Short_t **pdig= new Short_t*[nInputs]; //pointers to the expanded digits array
156 Int_t **ptr= new Int_t*[nInputs]; //pointers to teh expanded tracks array
157 Bool_t *active= new Bool_t[nInputs]; //flag for active input segments
160 //create digits array for given sectors
163 AliSimDigits ** digarr = new AliSimDigits*[nInputs];
164 for (Int_t i1=0;i1<nInputs; i1++)
168 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(i1));
169 gime = rl->GetLoader("TPCLoader");
170 gime->LoadSDigits("read");
171 TTree * treear = gime->TreeS();
175 cerr<<"AliTPCDigitizer: Input tree with SDigits not found in"
176 <<" input "<< i1<<endl;
180 if (treear->GetIndex()==0)
181 treear->BuildIndex("fSegmentID","fSegmentID");
182 treear->GetBranch("Segment")->SetAddress(&digarr[i1]);
186 //create branch's in TPC treeD
187 AliSimDigits * digrow = new AliSimDigits;
189 orl = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
190 ogime = orl->GetLoader("TPCLoader");
192 TTree * tree = ogime->TreeD();
195 ogime->MakeTree("D");
196 tree = ogime->TreeD();
198 tree->Branch("Segment","AliSimDigits",&digrow);
201 param->SetZeroSup(2);
203 Int_t zerosup = param->GetZeroSup();
204 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetDedxGainFactor();
205 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();
207 //Loop over segments of the TPC
209 for (Int_t segmentID=0; segmentID<param->GetNRowsTotal(); segmentID++)
212 if (!param->AdjustSectorRow(segmentID,sec,row))
214 cerr<<"AliTPC warning: invalid segment ID ! "<<segmentID<<endl;
217 AliTPCCalROC * gainROC = gainTPC->GetCalROC(sec); // pad gains per given sector
218 AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(sec); // noise per given sector
219 digrow->SetID(segmentID);
224 Bool_t digitize = kFALSE;
225 for (Int_t i=0;i<nInputs; i++)
228 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(i));
229 gime = rl->GetLoader("TPCLoader");
231 if (gime->TreeS()->GetEntryWithIndex(segmentID,segmentID) >= 0) {
232 digarr[i]->ExpandBuffer();
233 digarr[i]->ExpandTrackBuffer();
234 nrows = digarr[i]->GetNRows();
235 ncols = digarr[i]->GetNCols();
237 if (!fRegionOfInterest || (i == 0)) digitize = kTRUE;
241 if (fRegionOfInterest && !digitize) break;
243 if (!digitize) continue;
245 digrow->Allocate(nrows,ncols);
246 digrow->AllocateTrack(3);
249 Int_t label[1000]; //stack for 300 events
252 Int_t nElems = nrows*ncols;
254 for (Int_t i=0;i<nInputs; i++)
256 pdig[i] = digarr[i]->GetDigits();
257 ptr[i] = digarr[i]->GetTracks();
260 Short_t *pdig1= digrow->GetDigits();
261 Int_t *ptr1= digrow->GetTracks() ;
265 for (Int_t elem=0;elem<nElems; elem++)
271 for (Int_t i=0;i<nInputs; i++) if (active[i])
273 // q += digarr[i]->GetDigitFast(rows,col);
276 for (Int_t tr=0;tr<3;tr++)
278 // Int_t lab = digarr[i]->GetTrackIDFast(rows,col,tr);
279 Int_t lab = ptr[i][tr*nElems];
280 if ( (lab > 1) && *(pdig[i])>zerosup)
282 label[labptr]=lab+masks[i];
289 q/=16.; //conversion factor
290 Float_t gain = gainROC->GetValue(row,elem/nrows); // get gain for given - pad-row pad
292 //printf("problem\n");
295 Float_t noisePad = noiseROC->GetValue(row,elem/nrows);
296 // Float_t noise = gRandom->Gaus(0,param->GetNoise()*param->GetNoiseNormFac());
297 Float_t noise = pTPC->GetNoise();
302 if(q >= param->GetADCSat()) q = (Short_t)(param->GetADCSat() - 1);
303 //digrow->SetDigitFast((Short_t)q,rows,col);
305 for (Int_t tr=0;tr<3;tr++)
308 ptr1[tr*nElems] = label[tr];
315 digrow->CompresBuffer(1,zerosup);
316 digrow->CompresTrackBuffer(1);
318 if (fDebug>0) cerr<<sec<<"\t"<<row<<"\n";
319 } //for (Int_t n=0; n<param->GetNRowsTotal(); n++)
322 orl = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
323 ogime = orl->GetLoader("TPCLoader");
324 ogime->WriteDigits("OVERWRITE");
326 //fManager->GetTreeDTPC()->Write(0,TObject::kOverwrite);
329 for (Int_t i1=0;i1<nInputs; i1++) delete digarr[i1];
339 //------------------------------------------------------------------------
340 void AliTPCDigitizer::ExecSave(Option_t* option)
343 // merge input tree's with summable digits
344 //output stored in TreeTPCD
346 TString optionString = option;
347 if (!strcmp(optionString.Data(),"deb")) {
348 cout<<"AliTPCDigitizer::Exec: called with option deb "<<endl;
351 //get detector and geometry
352 AliRunLoader *rl, *orl;
353 AliLoader *gime, *ogime;
356 orl = AliRunLoader::GetRunLoader(fManager->GetOutputFolderName());
357 ogime = orl->GetLoader("TPCLoader");
359 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
360 gime = rl->GetLoader("TPCLoader");
363 AliRun* alirun = rl->GetAliRun();
365 AliTPC *pTPC = (AliTPC *) alirun->GetModule("TPC");
366 AliTPCParam * param = pTPC->GetParam();
367 pTPC->GenerNoise(500000); //create teble with noise
368 printf("noise %f \n", param->GetNoise()*param->GetNoiseNormFac());
370 Int_t nInputs = fManager->GetNinputs();
371 Int_t * masks = new Int_t[nInputs];
372 for (Int_t i=0; i<nInputs;i++)
373 masks[i]= fManager->GetMask(i);
375 AliSimDigits ** digarr = new AliSimDigits*[nInputs];
376 for (Int_t i1=0;i1<nInputs; i1++)
380 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(i1));
381 gime = rl->GetLoader("TPCLoader");
383 TTree * treear = gime->TreeS();
384 TBranch * br = treear->GetBranch("fSegmentID");
385 if (br) br->GetFile()->cd();
387 cerr<<" TPC - not existing input = \n"<<i1<<" ";
389 treear->GetBranch("Segment")->SetAddress(&digarr[i1]);
392 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
393 gime = rl->GetLoader("TPCLoader");
394 Stat_t nentries = gime->TreeS()->GetEntries();
397 //create branch's in TPC treeD
398 AliSimDigits * digrow = new AliSimDigits;
399 TTree * tree = ogime->TreeD();
401 tree->Branch("Segment","AliSimDigits",&digrow);
402 param->SetZeroSup(2);
404 Int_t zerosup = param->GetZeroSup();
405 //Loop over segments of the TPC
407 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetDedxGainFactor();
408 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();
409 for (Int_t n=0; n<nentries; n++) {
410 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(0));
411 gime = rl->GetLoader("TPCLoader");
412 gime->TreeS()->GetEvent(n);
414 digarr[0]->ExpandBuffer();
415 digarr[0]->ExpandTrackBuffer();
418 for (Int_t i=1;i<nInputs; i++){
419 // fManager->GetInputTreeTPCS(i)->GetEntryWithIndex(digarr[0]->GetID(),digarr[0]->GetID());
420 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(i));
421 gime = rl->GetLoader("TPCLoader");
422 gime->TreeS()->GetEntryWithIndex(digarr[0]->GetID(),digarr[0]->GetID());
423 digarr[i]->ExpandBuffer();
424 digarr[i]->ExpandTrackBuffer();
425 if ((digarr[0]->GetID()-digarr[i]->GetID())>0)
431 if (!param->AdjustSectorRow(digarr[0]->GetID(),sec,row)) {
432 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr[0]->GetID()<<endl;
436 AliTPCCalROC * gainROC = gainTPC->GetCalROC(sec); // pad gains per given sector
437 AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(sec); // noise per given sector
438 digrow->SetID(digarr[0]->GetID());
440 Int_t nrows = digarr[0]->GetNRows();
441 Int_t ncols = digarr[0]->GetNCols();
442 digrow->Allocate(nrows,ncols);
443 digrow->AllocateTrack(3);
446 Int_t label[1000]; //stack for 300 events
451 for (Int_t rows=0;rows<nrows; rows++){
452 for (Int_t col=0;col<ncols; col++){
457 for (Int_t i=0;i<nInputs; i++){
458 q += digarr[i]->GetDigitFast(rows,col);
461 for (Int_t tr=0;tr<3;tr++) {
462 Int_t lab = digarr[i]->GetTrackIDFast(rows,col,tr);
463 //Int_t lab = ptr[i][tr*nElems];
465 label[labptr]=lab+masks[i];
473 q/=16.; //conversion factor
474 // Float_t noise = gRandom->Gaus(0,param->GetNoise()*param->GetNoiseNormFac());
475 Float_t gain = gainROC->GetValue(row,col);
477 Float_t noisePad = noiseROC->GetValue(row, col);
479 Float_t noise = pTPC->GetNoise();
485 if(q >= param->GetADCSat()) q = (Short_t)(param->GetADCSat() - 1);
486 digrow->SetDigitFast((Short_t)q,rows,col);
487 // *pdig1 =Short_t(q);
488 for (Int_t tr=0;tr<3;tr++){
490 ((AliSimDigits*)digrow)->SetTrackIDFast(label[tr],rows,col,tr);
491 //ptr1[tr*nElems] = label[tr];
493 // ((AliSimDigits*)digrow)->SetTrackIDFast(-1,rows,col,tr);
494 // ptr1[tr*nElems] = 1;
502 digrow->CompresBuffer(1,zerosup);
503 digrow->CompresTrackBuffer(1);
505 if (fDebug>0) cerr<<sec<<"\t"<<row<<"\n";
507 // printf("end TPC merging - end -Tree %s\t%p\n",fManager->GetInputTreeH(0)->GetName(),fManager->GetInputTreeH(0)->GetListOfBranches()->At(3));
508 //fManager->GetTreeDTPC()->Write(0,TObject::kOverwrite);
509 ogime->WriteDigits("OVERWRITE");
511 for (Int_t i=1;i<nInputs; i++)
513 rl = AliRunLoader::GetRunLoader(fManager->GetInputFolderName(i));
514 gime = rl->GetLoader("TPCLoader");
515 gime->UnloadSDigits();
517 ogime->UnloadDigits();
520 for (Int_t i1=0;i1<nInputs; i1++) delete digarr[i1];