/************************************************************************** * Copyright(c) 1998-2000, ALICE Experiment at CERN, All rights reserved. * * * * Author: The ALICE Off-line Project. * * Contributors are mentioned in the code where appropriate. * * * * Permission to use, copy, modify and distribute this software and its * * documentation strictly for non-commercial purposes is hereby granted * * without fee, provided that the above copyright notice appears in all * * copies and that both the copyright notice and this permission notice * * appear in the supporting documentation. The authors make no claims * * about the suitability of this software for any purpose. It is * * provided "as is" without express or implied warranty. * **************************************************************************/ /* $Id$ */ /* Class for creating of the sumable digits and digits from MC data // The input : ideal signals (Hits->Diffusion->Attachment -Ideal signal) The output: "raw digits" Effect implemented: 1. Pad by pad gain map 2. Noise map 3. The dead channels identified - zerro noise for corresponding pads In this case the outpu equal zerro */ #include #include #include #include #include #include #include #include "AliTPCDigitizer.h" #include "AliTPC.h" #include "AliTPCParam.h" #include "AliTPCParamSR.h" #include "AliRun.h" #include "AliLoader.h" #include "AliPDG.h" #include "AliDigitizationInput.h" #include "AliSimDigits.h" #include "AliLog.h" #include "AliTPCcalibDB.h" #include "AliTPCCalPad.h" #include "AliTPCCalROC.h" #include "TTreeStream.h" #include "AliTPCReconstructor.h" #include using std::cout; using std::cerr; using std::endl; ClassImp(AliTPCDigitizer) //___________________________________________ AliTPCDigitizer::AliTPCDigitizer() :AliDigitizer(),fDebug(0), fDebugStreamer(0) { // // Default ctor - don't use it // } //___________________________________________ AliTPCDigitizer::AliTPCDigitizer(AliDigitizationInput* digInput) :AliDigitizer(digInput),fDebug(0), fDebugStreamer(0) { // // ctor which should be used // AliDebug(2,"(AliDigitizationInput* digInput) was processed"); if (AliTPCReconstructor::StreamLevel()>0) fDebugStreamer = new TTreeSRedirector("TPCDigitDebug.root","recreate"); } //------------------------------------------------------------------------ AliTPCDigitizer::~AliTPCDigitizer() { // Destructor if (fDebugStreamer) delete fDebugStreamer; } //------------------------------------------------------------------------ Bool_t AliTPCDigitizer::Init() { // Initialization return kTRUE; } //------------------------------------------------------------------------ void AliTPCDigitizer::Digitize(Option_t* option) { //DigitizeFast(option); DigitizeWithTailAndCrossTalk(option); } //------------------------------------------------------------------------ void AliTPCDigitizer::DigitizeFast(Option_t* option) { // merge input tree's with summable digits //output stored in TreeTPCD char s[100]; char ss[100]; TString optionString = option; if (!strcmp(optionString.Data(),"deb")) { cout<<"AliTPCDigitizer:::DigitizeFast called with option deb "<GetInputFolderName(0)); if (rl == 0x0) { Error("DigitizeFast","Can not find Run Loader for input 0. Can not proceed."); return; } rl->LoadgAlice(); rl->GetAliRun(); } AliTPC *pTPC = (AliTPC *) gAlice->GetModule("TPC"); AliTPCParam * param = pTPC->GetParam(); //sprintf(s,param->GetTitle()); snprintf(s,100,"%s",param->GetTitle()); //sprintf(ss,"75x40_100x60"); snprintf(ss,100,"75x40_100x60"); if(strcmp(s,ss)==0){ printf("2 pad-length geom hits with 3 pad-lenght geom digits...\n"); delete param; param=new AliTPCParamSR(); } else{ //sprintf(ss,"75x40_100x60_150x60"); snprintf(ss,100,"75x40_100x60_150x60"); if(strcmp(s,ss)!=0) { printf("No TPC parameters found...\n"); exit(2); } } pTPC->GenerNoise(500000); //create table with noise // Int_t nInputs = fDigInput->GetNinputs(); Int_t * masks = new Int_t[nInputs]; for (Int_t i=0; iGetMask(i); Short_t **pdig= new Short_t*[nInputs]; //pointers to the expanded digits array Int_t **ptr= new Int_t*[nInputs]; //pointers to the expanded tracks array Bool_t *active= new Bool_t[nInputs]; //flag for active input segments Char_t phname[100]; //create digits array for given sectors // make indexes // //create branch's in TPC treeD orl = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName()); ogime = orl->GetLoader("TPCLoader"); TTree * tree = ogime->TreeD(); AliSimDigits * digrow = new AliSimDigits; if (tree == 0x0) { ogime->MakeTree("D"); tree = ogime->TreeD(); } tree->Branch("Segment","AliSimDigits",&digrow); // AliSimDigits ** digarr = new AliSimDigits*[nInputs]; for (Int_t i1=0;i1GetInputFolderName(i1)); gime = rl->GetLoader("TPCLoader"); gime->LoadSDigits("read"); TTree * treear = gime->TreeS(); if (!treear) { cerr<<"AliTPCDigitizer: Input tree with SDigits not found in" <<" input "<< i1< *ph = (TParameter*)treear->GetUserInfo() ->FindObject("lhcphase0"); if(!ph){ cerr<<"AliTPCDigitizer: LHC phase not found in" <<" input "<< i1<GetUserInfo()->Add(new TParameter(phname,ph->GetVal())); // if (treear->GetIndex()==0) treear->BuildIndex("fSegmentID","fSegmentID"); treear->GetBranch("Segment")->SetAddress(&digarr[i1]); } // param->SetZeroSup(2); Int_t zerosup = param->GetZeroSup(); AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetDedxGainFactor(); AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise(); // //Loop over segments of the TPC for (Int_t segmentID=0; segmentIDGetNRowsTotal(); segmentID++) { Int_t sector, padRow; if (!param->AdjustSectorRow(segmentID,sector,padRow)) { cerr<<"AliTPC warning: invalid segment ID ! "<GetCalROC(sector); // pad gains per given sector AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(sector); // noise per given sector digrow->SetID(segmentID); Int_t nTimeBins = 0; Int_t nPads = 0; Bool_t digitize = kFALSE; for (Int_t i=0;iGetInputFolderName(i)); gime = rl->GetLoader("TPCLoader"); if (gime->TreeS()->GetEntryWithIndex(segmentID,segmentID) >= 0) { digarr[i]->ExpandBuffer(); digarr[i]->ExpandTrackBuffer(); nTimeBins = digarr[i]->GetNRows(); nPads = digarr[i]->GetNCols(); active[i] = kTRUE; if (!GetRegionOfInterest() || (i == 0)) digitize = kTRUE; } else { active[i] = kFALSE; } if (GetRegionOfInterest() && !digitize) break; } if (!digitize) continue; digrow->Allocate(nTimeBins,nPads); digrow->AllocateTrack(3); Float_t q=0; Int_t label[1000]; //stack for 300 events Int_t labptr = 0; Int_t nElems = nTimeBins*nPads; for (Int_t i=0;iGetDigits(); ptr[i] = digarr[i]->GetTracks(); } Short_t *pdig1= digrow->GetDigits(); Int_t *ptr1= digrow->GetTracks() ; for (Int_t elem=0;elemGetDigitFast(rows,col); q += *(pdig[i]); for (Int_t tr=0;tr<3;tr++) { // Int_t lab = digarr[i]->GetTrackIDFast(rows,col,tr); Int_t lab = ptr[i][tr*nElems]; if ( (lab > 1) && *(pdig[i])>zerosup) { label[labptr]=lab+masks[i]; labptr++; } } pdig[i]++; ptr[i]++; } q/=16.; //conversion factor Float_t gain = gainROC->GetValue(padRow,elem/nTimeBins); // get gain for given - pad-row pad //if (gain<0.5){ //printf("problem\n"); //} q*= gain; Float_t noisePad = noiseROC->GetValue(padRow,elem/nTimeBins); // Float_t noise = gRandom->Gaus(0,param->GetNoise()*param->GetNoiseNormFac()); Float_t noise = pTPC->GetNoise(); q+=noise*noisePad; q=TMath::Nint(q); if (q > zerosup) { if(q >= param->GetADCSat()) q = (Short_t)(param->GetADCSat() - 1); //digrow->SetDigitFast((Short_t)q,rows,col); *pdig1 =Short_t(q); for (Int_t tr=0;tr<3;tr++) { if (trGlitchFilter(); // digrow->CompresBuffer(1,zerosup); digrow->CompresTrackBuffer(1); tree->Fill(); if (fDebug>0) cerr<GetNRowsTotal(); n++) orl = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName()); ogime = orl->GetLoader("TPCLoader"); ogime->WriteDigits("OVERWRITE"); //fDigInput->GetTreeDTPC()->Write(0,TObject::kOverwrite); delete digrow; for (Int_t i1=0;i1GetOutputFolderName()); ogime = orl->GetLoader("TPCLoader"); rl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(0)); //gime = rl->GetLoader("TPCLoader"); rl->GetLoader("TPCLoader"); rl->LoadgAlice(); AliRun* alirun = rl->GetAliRun(); AliTPC *pTPC = (AliTPC *) alirun->GetModule("TPC"); AliTPCParam * param = pTPC->GetParam(); pTPC->GenerNoise(500000); //create teble with noise printf("noise %f \n", param->GetNoise()*param->GetNoiseNormFac()); // Int_t nInputs = fDigInput->GetNinputs(); // stupid protection... if (nInputs <= 0) return; // Int_t * masks = new Int_t[nInputs]; for (Int_t i=0; iGetMask(i); AliSimDigits ** digarr = new AliSimDigits*[nInputs]; for(Int_t ii=0;iiGetInputFolderName(i1)); gime = rl->GetLoader("TPCLoader"); TTree * treear = gime->TreeS(); // if (!treear) { cerr<<" TPC - not existing input = \n"<GetBranch("fSegmentID"); if (br) br->GetFile()->cd(); treear->GetBranch("Segment")->SetAddress(&digarr[i1]); } rl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(0)); gime = rl->GetLoader("TPCLoader"); Stat_t nentries = gime->TreeS()->GetEntries(); //create branch's in TPC treeD AliSimDigits * digrow = new AliSimDigits; TTree * tree = ogime->TreeD(); tree->Branch("Segment","AliSimDigits",&digrow); param->SetZeroSup(2); Int_t zerosup = param->GetZeroSup(); //Loop over segments of the TPC AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetDedxGainFactor(); AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise(); for (Int_t n=0; nGetInputFolderName(0)); gime = rl->GetLoader("TPCLoader"); gime->TreeS()->GetEvent(n); digarr[0]->ExpandBuffer(); digarr[0]->ExpandTrackBuffer(); for (Int_t i=1;iGetInputTreeTPCS(i)->GetEntryWithIndex(digarr[0]->GetID(),digarr[0]->GetID()); rl = AliRunLoader::GetRunLoader(fDigInput->GetInputFolderName(i)); gime = rl->GetLoader("TPCLoader"); gime->TreeS()->GetEntryWithIndex(digarr[0]->GetID(),digarr[0]->GetID()); digarr[i]->ExpandBuffer(); digarr[i]->ExpandTrackBuffer(); if ((digarr[0]->GetID()-digarr[i]->GetID())>0) printf("problem\n"); } Int_t sector, padRow; if (!param->AdjustSectorRow(digarr[0]->GetID(),sector,padRow)) { cerr<<"AliTPC warning: invalid segment ID ! "<GetID()<GetCalROC(sector); // pad gains per given sector AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(sector); // noise per given sector digrow->SetID(digarr[0]->GetID()); Int_t nTimeBins = digarr[0]->GetNRows(); Int_t nPads = digarr[0]->GetNCols(); digrow->Allocate(nTimeBins,nPads); digrow->AllocateTrack(3); Float_t q=0; Int_t label[1000]; //stack for 300 events Int_t labptr = 0; for (Int_t iTimeBin=0;iTimeBinGetDigitFast(iTimeBin,iPad); //q += *(pdig[i]); for (Int_t tr=0;tr<3;tr++) { Int_t lab = digarr[i]->GetTrackIDFast(iTimeBin,iPad,tr); //Int_t lab = ptr[i][tr*nElems]; if ( (lab > 1) ) { label[labptr]=lab+masks[i]; labptr++; } } // pdig[i]++; //ptr[i]++; } q/=16.; //conversion factor // Float_t noise = gRandom->Gaus(0,param->GetNoise()*param->GetNoiseNormFac()); Float_t gain = gainROC->GetValue(padRow,iPad); q*= gain; Float_t noisePad = noiseROC->GetValue(padRow, iPad); Float_t noise = pTPC->GetNoise(); q+=noise*noisePad; // // here we can get digits from past and add signal // // //for (Int_t jTimeBin=0; jTimeBin zerosup){ if(q >= param->GetADCSat()) q = (Short_t)(param->GetADCSat() - 1); digrow->SetDigitFast((Short_t)q,iTimeBin,iPad); // *pdig1 =Short_t(q); for (Int_t tr=0;tr<3;tr++){ if (trSetTrackIDFast(label[tr],iTimeBin,iPad,tr); //ptr1[tr*nElems] = label[tr]; //else // ((AliSimDigits*)digrow)->SetTrackIDFast(-1,iTimeBin,iPad,tr); // ptr1[tr*nElems] = 1; } } //pdig1++; //ptr1++; } } digrow->CompresBuffer(1,zerosup); digrow->CompresTrackBuffer(1); tree->Fill(); if (fDebug>0) cerr<GetInputTreeH(0)->GetName(),fDigInput->GetInputTreeH(0)->GetListOfBranches()->At(3)); //fDigInput->GetTreeDTPC()->Write(0,TObject::kOverwrite); ogime->WriteDigits("OVERWRITE"); for (Int_t i=1;iGetInputFolderName(i)); gime = rl->GetLoader("TPCLoader"); gime->UnloadSDigits(); } ogime->UnloadDigits(); delete digrow; for (Int_t i1=0;i1GetRecoParam(0); AliDebug(1, Form(" recoParam->GetCrosstalkCorrection() = %f", recoParam->GetCrosstalkCorrection())); AliDebug(1,Form(" recoParam->GetUseIonTailCorrection() = %d ", recoParam->GetUseIonTailCorrection())); Int_t nROCs = 72; char s[100]; char ss[100]; TString optionString = option; if (!strcmp(optionString.Data(),"deb")) { cout<<"AliTPCDigitizer:::DigitizeFast called with option deb "<GetInputFolderName(0)); if (rl == 0x0) { Error("DigitizeFast","Can not find Run Loader for input 0. Can not proceed."); return; } rl->LoadgAlice(); rl->GetAliRun(); } AliTPC *pTPC = (AliTPC *) gAlice->GetModule("TPC"); AliTPCParam * param = pTPC->GetParam(); //sprintf(s,param->GetTitle()); snprintf(s,100,"%s",param->GetTitle()); //sprintf(ss,"75x40_100x60"); snprintf(ss,100,"75x40_100x60"); if(strcmp(s,ss)==0){ printf("2 pad-length geom hits with 3 pad-lenght geom digits...\n"); delete param; param=new AliTPCParamSR(); } else { //sprintf(ss,"75x40_100x60_150x60"); snprintf(ss,100,"75x40_100x60_150x60"); if(strcmp(s,ss)!=0) { printf("No TPC parameters found...\n"); exit(2); } } pTPC->GenerNoise(500000); //create table with noise // Int_t nInputs = fDigInput->GetNinputs(); Int_t * masks = new Int_t[nInputs]; for (Int_t i=0; iGetMask(i); Short_t **pdig= new Short_t*[nInputs]; //pointers to the expanded digits array Int_t **ptr= new Int_t*[nInputs]; //pointers to the expanded tracks array Bool_t *active= new Bool_t[nInputs]; //flag for active input segments Char_t phname[100]; //create digits array for given sectors // make indexes // //create branch's in TPC treeD orl = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName()); ogime = orl->GetLoader("TPCLoader"); TTree * tree = ogime->TreeD(); AliSimDigits * digrow = new AliSimDigits; if (tree == 0x0) { ogime->MakeTree("D"); tree = ogime->TreeD(); } tree->Branch("Segment","AliSimDigits",&digrow); // AliSimDigits ** digarr = new AliSimDigits*[nInputs]; for (Int_t i1=0;i1GetInputFolderName(i1)); gime = rl->GetLoader("TPCLoader"); gime->LoadSDigits("read"); TTree * treear = gime->TreeS(); if (!treear) { cerr<<"AliTPCDigitizer: Input tree with SDigits not found in" <<" input "<< i1< *ph = (TParameter*)treear->GetUserInfo() ->FindObject("lhcphase0"); if(!ph) { cerr<<"AliTPCDigitizer: LHC phase not found in" <<" input "<< i1<GetUserInfo()->Add(new TParameter(phname,ph->GetVal())); // if (treear->GetIndex()==0) treear->BuildIndex("fSegmentID","fSegmentID"); treear->GetBranch("Segment")->SetAddress(&digarr[i1]); } // // zero supp, take gain and noise map of TPC from OCDB param->SetZeroSup(2); Int_t zerosup = param->GetZeroSup(); AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetDedxGainFactor(); AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise(); // // Cache the ion tail objects form OCDB // TObjArray *ionTailArr = (TObjArray*)AliTPCcalibDB::Instance()->GetIonTailArray(); if (!ionTailArr) {AliFatal("TPC - Missing IonTail OCDB object");} // TObject *rocFactorIROC = ionTailArr->FindObject("factorIROC"); // TObject *rocFactorOROC = ionTailArr->FindObject("factorOROC"); // Float_t factorIROC = (atof(rocFactorIROC->GetTitle())); // Float_t factorOROC = (atof(rocFactorOROC->GetTitle())); Int_t nIonTailBins =0; TObjArray timeResFunc(nROCs); for (Int_t isec = 0;isecGetTailcancelationGraphs(isec,graphRes,trfIndexArr)) continue; // fill all TGraphErrors of trfs (time response function) of a given sector to a TObjArray TObjArray *timeResArr = new TObjArray(20); timeResArr -> SetOwner(kTRUE); for (Int_t ires = 0;ires<20;ires++) timeResArr->AddAt(graphRes[ires],ires); timeResFunc.AddAt(timeResArr,isec); // Fill all trfs into a single TObjArray nIonTailBins = graphRes[3]->GetN(); } // // 1.) Make first loop to calculate mean amplitude per pad per segment for cross talk // TObjArray crossTalkSignalArray(nROCs); // for 36 sectors TVectorD * qTotSectorOld = new TVectorD(nROCs); Float_t qTotTPC = 0.; Float_t qTotPerSector = 0.; Int_t nTimeBinsAll = 1100; Int_t nWireSegments=11; // 1.a) crorstalk matrix initialization for (Int_t sector=0; sectorGetNRowsTotal(); globalRowID++) { Int_t sector, padRow; if (!param->AdjustSectorRow(globalRowID,sector,padRow)) { cerr<<"AliTPC warning: invalid segment ID ! "<GetWireSegment(sector,padRow); Float_t nPadsPerSegment = (Float_t)(param->GetNPadsPerSegment(wireSegmentID)); // structure with mean signal per pad to be filled for each timebin in first loop (11 anodeWireSegment and 1100 timebin) TMatrixD &crossTalkSignal = *((TMatrixD*)crossTalkSignalArray.At(sector)); AliTPCCalROC * gainROC = gainTPC->GetCalROC(sector); // pad gains per given sector // digrow->SetID(globalRowID); Int_t nTimeBins = 0; Int_t nPads = 0; Bool_t digitize = kFALSE; for (Int_t i=0;iGetInputFolderName(i)); gime = rl->GetLoader("TPCLoader"); if (gime->TreeS()->GetEntryWithIndex(globalRowID,globalRowID) >= 0) { digarr[i]->ExpandBuffer(); digarr[i]->ExpandTrackBuffer(); nTimeBins = digarr[i]->GetNRows(); nPads = digarr[i]->GetNCols(); active[i] = kTRUE; if (!GetRegionOfInterest() || (i == 0)) digitize = kTRUE; } else { active[i] = kFALSE; } if (GetRegionOfInterest() && !digitize) break; } if (!digitize) continue; //digrow->Allocate(nTimeBins,nPads); Float_t q = 0; Int_t labptr = 0; Int_t nElems = nTimeBins*nPads; // element is a unit of a given row's pad-timebin space for (Int_t i=0;iGetDigits(); } // // loop over elements i.e pad-timebin space of a row, "padNumber=elem/nTimeBins", "timeBin=elem%nTimeBins" // for (Int_t elem=0;elemGetValue(padRow,padNumber); // get gain for given - pad-row pad q*= gain; crossTalkSignal[wireSegmentID][timeBin]+= q/nPadsPerSegment; // Qtot per segment for a given timebin qTotSectorOld -> GetMatrixArray()[sector] += q; // Qtot for each sector qTotTPC += q; // Qtot for whole TPC } // end of q loop } // end of global row loop // // 1.b) Dump the content of the crossTalk signal to the debug stremer - to be corealted later with the same crosstalk correction // assumed during reconstruction // if (AliTPCReconstructor::StreamLevel()==1) { // // dump the crosstalk matrices to tree for further investigation // a.) to estimate fluctuation of pedestal in indiviula wire segments // b.) to check correlation between regions // c.) to check relative conribution of signal below threshold to crosstalk for (Int_t isector=0; isectorAt(isector+nROCs); TVectorD vecAll(crossTalkMatrix->GetNrows()); //TVectorD vecBelow(crossTalkMatrix->GetNrows()); // for (Int_t itime=0; itimeGetNcols(); itime++){ for (Int_t iwire=0; iwireGetNrows(); iwire++){ vecAll[iwire]=(*crossTalkMatrix)(iwire,itime); //vecBelow[iwire]=(*crossTalkMatrixBelow)(iwire,itime); } (*fDebugStreamer)<<"crosstalkMatrix"<< "sector="<GetNRowsTotal(); globalRowID++) { Int_t sector, padRow; if (!param->AdjustSectorRow(globalRowID,sector,padRow)) { cerr<<"AliTPC warning: invalid segment ID ! "<At(1); Int_t wireSegmentID = param->GetWireSegment(sector,padRow); Float_t nPadsPerSegment = (Float_t)(param->GetNPadsPerSegment(wireSegmentID)); // const Float_t ampfactor = (sector<36)?factorIROC:factorOROC; // factor for the iontail which is ROC type dependent AliTPCCalROC * gainROC = gainTPC->GetCalROC(sector); // pad gains per given sector AliTPCCalROC * noiseROC = noiseTPC->GetCalROC(sector); // noise per given sector digrow->SetID(globalRowID); Int_t nTimeBins = 0; Int_t nPads = 0; Bool_t digitize = kFALSE; for (Int_t i=0;iGetInputFolderName(i)); gime = rl->GetLoader("TPCLoader"); if (gime->TreeS()->GetEntryWithIndex(globalRowID,globalRowID) >= 0) { digarr[i]->ExpandBuffer(); digarr[i]->ExpandTrackBuffer(); nTimeBins = digarr[i]->GetNRows(); nPads = digarr[i]->GetNCols(); active[i] = kTRUE; if (!GetRegionOfInterest() || (i == 0)) digitize = kTRUE; } else { active[i] = kFALSE; } if (GetRegionOfInterest() && !digitize) break; } if (!digitize) continue; digrow->Allocate(nTimeBins,nPads); digrow->AllocateTrack(3); Int_t localPad = 0; Float_t q = 0.; Float_t qXtalk = 0.; Float_t qIonTail = 0.; Float_t qOrig = 0.; Int_t label[1000]; //stack for 300 events Int_t labptr = 0; Int_t nElems = nTimeBins*nPads; // element is a unit of a given row's pad-timebin space for (Int_t i=0;iGetDigits(); ptr[i] = digarr[i]->GetTracks(); } Short_t *pdig1= digrow->GetDigits(); Int_t *ptr1= digrow->GetTracks() ; // loop over elements i.e pad-timebin space of a row for (Int_t elem=0;elem 1) && *(pdig[i])>zerosup) { label[labptr]=lab+masks[i]; labptr++; } } pdig[i]++; ptr[i]++; } Int_t padNumber = elem/nTimeBins; Int_t timeBin = elem%nTimeBins; localPad = padNumber-nPads/2; Float_t gain = gainROC->GetValue(padRow,padNumber); // get gain for given - pad-row pad //if (gain<0.5){ //printf("problem\n"); //} q*= gain; qOrig = q; Float_t noisePad = noiseROC->GetValue(padRow,padNumber); Float_t noise = pTPC->GetNoise()*noisePad; if ( (q/16.+noise)> zerosup){ // Crosstalk correction qXtalk = (*(TMatrixD*)crossTalkSignalArray.At(sector))[wireSegmentID][timeBin]; qTotPerSector = qTotSectorOld -> GetMatrixArray()[sector]; // Ion tail correction: being elem=padNumber*nTimeBins+timeBin; Int_t lowerElem=elem-nIonTailBins; Int_t zeroElem =(elem/nTimeBins)*nTimeBins; if (zeroElem<0) zeroElem=0; if (lowerElem0 && recoParam->GetUseIonTailCorrection()){ for (Int_t i=0;iGetDigits(); if (padNumber==0) continue; if (padNumber>=nPads-1) continue; for (Int_t dpad=-1; dpad<=1; dpad++){ // calculate iontail due signals from neigborhood pads for (Int_t celem=elem-1; celem>lowerElem; celem--){ Int_t celemPad=celem+dpad*nTimeBins; Double_t qCElem=pdigC[celemPad]; if ( qCElem<=0) continue; //here we substract ion tail Double_t COG=0; if (celemPad-nTimeBins>nTimeBins && celemPad+nTimeBinsAt(indexTRFPRF); if (graphTRFPRF==NULL) continue; // here we should get index and point of TRF corresponding to given COG position if (graphTRFPRF->GetY()[elem-celem]<0)qIonTail+=qCElem*graphTRFPRF->GetY()[elem-celem]; } } } } } // q -= qXtalk*recoParam->GetCrosstalkCorrection(); q+=qIonTail; q/=16.; //conversion factor q+=noise; q=TMath::Nint(q); // round to the nearest integer // fill info for degugging if (AliTPCReconstructor::StreamLevel()==1 && qOrig > zerosup ) { TTreeSRedirector &cstream = *fDebugStreamer; UInt_t uid = AliTPCROC::GetTPCUniqueID(sector, padRow, padNumber,timeBin); cstream <<"ionTailXtalk"<< "uid="< zerosup){ if(q >= param->GetADCSat()) q = (Short_t)(param->GetADCSat() - 1); //digrow->SetDigitFast((Short_t)q,rows,col); *pdig1 =Short_t(q); for (Int_t tr=0;tr<3;tr++) { if (tr zerosup ) { cout << " sector = " << sector << " row = " << padRow << " localPad = " << localPad << " SegmentID = " << wireSegmentID << " nPadsPerSegment = " << nPadsPerSegment << endl; cout << " qXtalk = " << qXtalk ; cout << " qOrig = " << qOrig ; cout << " q = " << q ; cout << " qsec = " << qTotPerSector ; cout << " qTotTPC "<< qTotTPC << endl; } // // glitch filter // digrow->GlitchFilter(); // digrow->CompresBuffer(1,zerosup); digrow->CompresTrackBuffer(1); tree->Fill(); if (fDebug>0) cerr<GetNRowsTotal(); n++) orl = AliRunLoader::GetRunLoader(fDigInput->GetOutputFolderName()); ogime = orl->GetLoader("TPCLoader"); ogime->WriteDigits("OVERWRITE"); //fDigInput->GetTreeDTPC()->Write(0,TObject::kOverwrite); delete digrow; for (Int_t i1=0;i1