]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCclustererKr.cxx
move minuit initialization to the unfolding method since it leaks and it is not used...
[u/mrichter/AliRoot.git] / TPC / AliTPCclustererKr.cxx
CommitLineData
e2a4d72c 1/**************************************************************************\r
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *\r
3 * *\r
4 * Author: The ALICE Off-line Project. *\r
5 * Contributors are mentioned in the code where appropriate. *\r
6 * *\r
7 * Permission to use, copy, modify and distribute this software and its *\r
8 * documentation strictly for non-commercial purposes is hereby granted *\r
9 * without fee, provided that the above copyright notice appears in all *\r
10 * copies and that both the copyright notice and this permission notice *\r
11 * appear in the supporting documentation. The authors make no claims *\r
12 * about the suitability of this software for any purpose. It is *\r
13 * provided "as is" without express or implied warranty. *\r
14 **************************************************************************/\r
15\r
16/* $Id: AliTPCclustererKr.cxx,v 1.7 2008/02/06 17:24:53 matyja Exp $ */\r
17\r
18//-----------------------------------------------------------------\r
19// Implementation of the TPC Kr cluster class\r
20//\r
21// Origin: Adam Matyja, INP PAN, adam.matyja@ifj.edu.pl\r
22//-----------------------------------------------------------------\r
23\r
24/*\r
25Instruction - how to use that:\r
26There are two macros prepared. One is for preparing clusters from MC \r
27samples: FindKrClusters.C. The output is kept in TPC.RecPoints.root.\r
28The other macro is prepared for data analysis: FindKrClustersRaw.C. \r
29The output is created for each processed file in root file named adc.root. \r
30For each data subsample the same named file is created. So be careful \r
31do not overwrite them. \r
32\r
33Additional selection criteria to select the GOLD cluster\r
34Example:\r
35// open file with clusters\r
36TFile f("Krypton.root");\r
37TTree * tree = (TTree*)f.Get("Kr")\r
38TCut cutR0("cutR0","fADCcluster/fSize<100"); // adjust it according v seetings - \r
39TCut cutR1("cutR1","fADCcluster/fSize>7"); // cosmic tracks and noise removal\r
40TCut cutR2("cutR2","fMax.fAdc/fADCcluster<0.2"); // digital noise removal\r
41TCut cutR3("cutR3","fMax.fAdc/fADCcluster>0.01"); // noise removal\r
42TCut cutS1("cutS1","fSize<200"); // adjust it according v seetings - cosmic tracks\r
43TCut cutAll = cutR0+cutR1+cutR2+cutR3+cutS1;\r
44This values are typical values to be applied in selectors\r
45\r
46\r
47*\r
48**** MC ****\r
49*\r
50\r
51To run clusterizaton for MC type:\r
52.x FindKrClusters.C\r
53\r
54If you don't want to use the standard selection criteria then you \r
55have to do following:\r
56\r
57// load the standard setup\r
58AliRunLoader* rl = AliRunLoader::Open("galice.root");\r
59AliTPCLoader *tpcl = (AliTPCLoader*)rl->GetLoader("TPCLoader");\r
60tpcl->LoadDigits();\r
61rl->LoadgAlice();\r
62gAlice=rl->GetAliRun();\r
63TDirectory *cwd = gDirectory;\r
64AliTPCv4 *tpc = (AliTPCv4*)gAlice->GetDetector("TPC");\r
65Int_t ver = tpc->IsVersion();\r
66rl->CdGAFile();\r
67AliTPCParam *param=(AliTPCParamSR *)gDirectory->Get("75x40_100x60_150x60");\r
68AliTPCDigitsArray *digarr=new AliTPCDigitsArray;\r
69digarr->Setup(param);\r
70cwd->cd();\r
71\r
72//loop over events\r
73Int_t nevmax=rl->GetNumberOfEvents();\r
74for(Int_t nev=0;nev<nevmax ;nev++){\r
75 rl->GetEvent(nev);\r
76 TTree* input_tree= tpcl->TreeD();//load tree with digits\r
77 digarr->ConnectTree(input_tree);\r
78 TTree *output_tree =tpcl->TreeR();//load output tree\r
79\r
80 AliTPCclustererKr *clusters = new AliTPCclustererKr();\r
81 clusters->SetParam(param);\r
82 clusters->SetInput(input_tree);\r
83 clusters->SetOutput(output_tree);\r
84 clusters->SetDigArr(digarr);\r
85 \r
86//If you want to change the cluster finder parameters for MC there are \r
87//several of them:\r
88\r
89//1. signal threshold (everything below the given number is treated as 0)\r
90 clusters->SetMinAdc(3);\r
91\r
92//2. number of neighbouring timebins to be considered\r
93 clusters->SetMinTimeBins(2);\r
94\r
95//3. distance of the cluster center to the center of a pad in pad-padrow plane \r
96//(in cm). Remenber that this is still quantified by pad size.\r
97 clusters->SetMaxPadRangeCm(2.5);\r
98\r
99//4. distance of the cluster center to the center of a padrow in pad-padrow \r
100//plane (in cm). Remenber that this is still quantified by pad size.\r
101 clusters->SetMaxRowRangeCm(3.5);\r
102\r
103//5. distance of the cluster center to the max time bin on a pad (in tackts)\r
104//ie. fabs(centerT - time)<7\r
105 clusters->SetMaxTimeRange(7);\r
106\r
107//6. cut reduce peak at 0. There are noises which appear mostly as two \r
108//timebins on one pad.\r
109 clusters->SetValueToSize(3.5);\r
110\r
111\r
112 clusters->finderIO();\r
113 tpcl->WriteRecPoints("OVERWRITE");\r
114}\r
115delete rl;//cleans everything\r
116\r
117*\r
118********* DATA *********\r
119*\r
120\r
121To run clusterizaton for DATA for file named raw_data.root type:\r
122.x FindKrClustersRaw.C("raw_data.root")\r
123\r
124If you want to change some criteria do the following:\r
125\r
126//\r
127// remove Altro warnings\r
128//\r
129AliLog::SetClassDebugLevel("AliTPCRawStream",-5);\r
130AliLog::SetClassDebugLevel("AliRawReaderDate",-5);\r
131AliLog::SetClassDebugLevel("AliTPCAltroMapping",-5);\r
132AliLog::SetModuleDebugLevel("RAW",-5);\r
133\r
134//\r
135// Get database with noises\r
136//\r
137// char *ocdbpath = gSystem->Getenv("OCDB_PATH");\r
138char *ocdbpath ="local:///afs/cern.ch/alice/tpctest/OCDB";\r
139if (ocdbpath==0){\r
140ocdbpath="alien://folder=/alice/data/2007/LHC07w/OCDB/";\r
141}\r
142AliCDBManager * man = AliCDBManager::Instance();\r
143man->SetDefaultStorage(ocdbpath);\r
144man->SetRun(0);\r
145AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();\r
146AliTPCAltroMapping** mapping =AliTPCcalibDB::Instance()->GetMapping();\r
147\r
148//define tree\r
149TFile *hfile=new TFile("adc.root","RECREATE","ADC file");\r
150// Create a ROOT Tree\r
151TTree *mytree = new TTree("Kr","Krypton cluster tree");\r
152\r
153//define infput file\r
154const char *fileName="data.root";\r
155AliRawReader *reader = new AliRawReaderRoot(fileName);\r
156//AliRawReader *reader = new AliRawReaderDate(fileName);\r
157reader->Reset();\r
158AliAltroRawStreamFast* stream = new AliAltroRawStreamFast(reader);\r
159stream->SelectRawData("TPC");\r
160\r
161//one general output\r
162AliTPCclustererKr *clusters = new AliTPCclustererKr();\r
163clusters->SetOutput(mytree);\r
164clusters->SetRecoParam(0);//standard reco parameters\r
165AliTPCParamSR *param=new AliTPCParamSR();\r
166clusters->SetParam(param);//TPC parameters(sectors, timebins, etc.)\r
167\r
168//set cluster finder parameters (from data):\r
169//1. zero suppression parameter\r
170 clusters->SetZeroSup(param->GetZeroSup());\r
171\r
172//2. first bin\r
173 clusters->SetFirstBin(60);\r
174\r
175//3. last bin\r
176 clusters->SetLastBin(950);\r
177\r
178//4. maximal noise\r
179 clusters->SetMaxNoiseAbs(2);\r
180\r
181//5. maximal amount of sigma of noise\r
182 clusters->SetMaxNoiseSigma(3);\r
183\r
184//The remaining parameters are the same paramters as for MC (see MC section \r
185//points 1-6)\r
186 clusters->SetMinAdc(3);\r
187 clusters->SetMinTimeBins(2);\r
188 clusters->SetMaxPadRangeCm(2.5);\r
189 clusters->SetMaxRowRangeCm(3.5);\r
190 clusters->SetMaxTimeRange(7);\r
191 clusters->SetValueToSize(3.5);\r
192\r
193while (reader->NextEvent()) {\r
194 clusters->FinderIO(reader);\r
195}\r
196\r
197hfile->Write();\r
198hfile->Close();\r
199delete stream;\r
200\r
201\r
202*/\r
203\r
204#include "AliTPCclustererKr.h"\r
205#include "AliTPCclusterKr.h"\r
206//#include <vector>\r
207#include <list>\r
208#include "TObject.h"\r
209#include "AliPadMax.h"\r
210#include "AliSimDigits.h"\r
211#include "AliTPCv4.h"\r
212#include "AliTPCParam.h"\r
213#include "AliTPCDigitsArray.h"\r
214#include "AliTPCvtpr.h"\r
215#include "AliTPCClustersRow.h"\r
216#include "TTree.h"\r
217#include "TH1F.h"\r
218#include "TH2F.h"\r
219#include "TTreeStream.h"\r
220\r
221#include "AliTPCTransform.h"\r
222\r
223//used in raw data finder\r
224#include "AliTPCROC.h"\r
225#include "AliTPCCalPad.h"\r
226#include "AliTPCAltroMapping.h"\r
227#include "AliTPCcalibDB.h"\r
228#include "AliTPCRawStream.h"\r
98b76430 229#include "AliTPCRawStreamV3.h"\r
e2a4d72c 230#include "AliTPCRecoParam.h"\r
231#include "AliTPCReconstructor.h"\r
232#include "AliRawReader.h"\r
233#include "AliTPCCalROC.h"\r
98b76430 234#include "AliRawEventHeaderBase.h"\r
e2a4d72c 235\r
236ClassImp(AliTPCclustererKr)\r
237\r
238\r
239AliTPCclustererKr::AliTPCclustererKr()\r
240 :TObject(),\r
241 fRawData(kFALSE),\r
242 fInput(0),\r
243 fOutput(0),\r
244 fParam(0),\r
245 fDigarr(0),\r
246 fRecoParam(0),\r
247 fZeroSup(2),\r
248 fFirstBin(60),\r
249 fLastBin(950),\r
250 fMaxNoiseAbs(2),\r
251 fMaxNoiseSigma(3),\r
252 fMinAdc(3),\r
253 fMinTimeBins(2),\r
254// fMaxPadRange(4),\r
255// fMaxRowRange(3),\r
256 fMaxTimeRange(7),\r
257 fValueToSize(3.5),\r
258 fMaxPadRangeCm(2.5),\r
259 fMaxRowRangeCm(3.5),\r
260 fIsolCut(3),\r
261 fDebugLevel(-1),\r
262 fHistoRow(0),\r
263 fHistoPad(0),\r
264 fHistoTime(0),\r
f5178f60 265 fHistoRowPad(0),\r
266 fTimeStamp(0),\r
267 fRun(0)\r
e2a4d72c 268{\r
269//\r
270// default constructor\r
271//\r
272}\r
273\r
274AliTPCclustererKr::AliTPCclustererKr(const AliTPCclustererKr &param)\r
275 :TObject(),\r
276 fRawData(kFALSE),\r
277 fInput(0),\r
278 fOutput(0),\r
279 fParam(0),\r
280 fDigarr(0),\r
281 fRecoParam(0),\r
282 fZeroSup(2),\r
283 fFirstBin(60),\r
284 fLastBin(950),\r
285 fMaxNoiseAbs(2),\r
286 fMaxNoiseSigma(3),\r
287 fMinAdc(3),\r
288 fMinTimeBins(2),\r
289// fMaxPadRange(4),\r
290// fMaxRowRange(3),\r
291 fMaxTimeRange(7),\r
292 fValueToSize(3.5),\r
293 fMaxPadRangeCm(2.5),\r
294 fMaxRowRangeCm(3.5),\r
295 fIsolCut(3),\r
296 fDebugLevel(-1),\r
297 fHistoRow(0),\r
298 fHistoPad(0),\r
299 fHistoTime(0),\r
98b76430 300 fHistoRowPad(0),\r
f5178f60 301 fTimeStamp(0),\r
302 fRun(0)\r
e2a4d72c 303{\r
304//\r
305// copy constructor\r
306//\r
307 fParam = param.fParam;\r
308 fRecoParam = param.fRecoParam;\r
309 fRawData = param.fRawData;\r
310 fInput = param.fInput ;\r
311 fOutput = param.fOutput;\r
312 fDigarr = param.fDigarr;\r
313 fZeroSup = param.fZeroSup ;\r
314 fFirstBin = param.fFirstBin ;\r
315 fLastBin = param.fLastBin ;\r
316 fMaxNoiseAbs = param.fMaxNoiseAbs ;\r
317 fMaxNoiseSigma = param.fMaxNoiseSigma ;\r
318 fMinAdc = param.fMinAdc;\r
319 fMinTimeBins = param.fMinTimeBins;\r
320// fMaxPadRange = param.fMaxPadRange ;\r
321// fMaxRowRange = param.fMaxRowRange ;\r
322 fMaxTimeRange = param.fMaxTimeRange;\r
323 fValueToSize = param.fValueToSize;\r
324 fMaxPadRangeCm = param.fMaxPadRangeCm;\r
325 fMaxRowRangeCm = param.fMaxRowRangeCm;\r
326 fIsolCut = param.fIsolCut;\r
327 fDebugLevel = param.fDebugLevel;\r
328 fHistoRow = param.fHistoRow ;\r
329 fHistoPad = param.fHistoPad ;\r
330 fHistoTime = param.fHistoTime;\r
331 fHistoRowPad = param.fHistoRowPad;\r
98b76430 332 fTimeStamp = param.fTimeStamp;\r
f5178f60 333 fRun = param.fRun;\r
e2a4d72c 334\r
335} \r
336\r
337AliTPCclustererKr & AliTPCclustererKr::operator = (const AliTPCclustererKr & param)\r
338{\r
339 //\r
340 // assignment operator\r
341 //\r
342 fParam = param.fParam;\r
343 fRecoParam = param.fRecoParam;\r
344 fRawData = param.fRawData;\r
345 fInput = param.fInput ;\r
346 fOutput = param.fOutput;\r
347 fDigarr = param.fDigarr;\r
348 fZeroSup = param.fZeroSup ;\r
349 fFirstBin = param.fFirstBin ;\r
350 fLastBin = param.fLastBin ;\r
351 fMaxNoiseAbs = param.fMaxNoiseAbs ;\r
352 fMaxNoiseSigma = param.fMaxNoiseSigma ;\r
353 fMinAdc = param.fMinAdc;\r
354 fMinTimeBins = param.fMinTimeBins;\r
355// fMaxPadRange = param.fMaxPadRange ;\r
356// fMaxRowRange = param.fMaxRowRange ;\r
357 fMaxTimeRange = param.fMaxTimeRange;\r
358 fValueToSize = param.fValueToSize;\r
359 fMaxPadRangeCm = param.fMaxPadRangeCm;\r
360 fMaxRowRangeCm = param.fMaxRowRangeCm;\r
361 fIsolCut = param.fIsolCut;\r
362 fDebugLevel = param.fDebugLevel;\r
363 fHistoRow = param.fHistoRow ;\r
364 fHistoPad = param.fHistoPad ;\r
365 fHistoTime = param.fHistoTime;\r
366 fHistoRowPad = param.fHistoRowPad;\r
98b76430 367 fTimeStamp = param.fTimeStamp;\r
f5178f60 368 fRun = param.fRun;\r
e2a4d72c 369 return (*this);\r
370}\r
371\r
372AliTPCclustererKr::~AliTPCclustererKr()\r
373{\r
374 //\r
375 // destructor\r
376 //\r
377 delete fOutput;\r
378}\r
379\r
380void AliTPCclustererKr::SetRecoParam(AliTPCRecoParam *recoParam)\r
381{\r
382 //\r
383 // set reconstruction parameters\r
384 //\r
385 if (recoParam) {\r
386 fRecoParam = recoParam;\r
387 }else{\r
388 //set default parameters if not specified\r
389 fRecoParam = AliTPCReconstructor::GetRecoParam();\r
390 if (!fRecoParam) fRecoParam = AliTPCRecoParam::GetLowFluxParam();\r
391 }\r
392 return;\r
393}\r
394\r
395\r
396////____________________________________________________________________________\r
397//// I/O\r
398void AliTPCclustererKr::SetInput(TTree * tree)\r
399{\r
400 //\r
401 // set input tree with digits\r
402 //\r
403 fInput = tree; \r
404 if (!fInput->GetBranch("Segment")){\r
405 cerr<<"AliTPCclusterKr::FindClusterKr(): no proper input tree !\n";\r
406 fInput=0;\r
407 return;\r
408 }\r
409}\r
410\r
411void AliTPCclustererKr::SetOutput(TTree * /*tree*/) \r
412{\r
413 //\r
414 //dummy\r
415 //\r
416 fOutput = new TTreeSRedirector("Krypton.root");\r
417}\r
418\r
419////____________________________________________________________________________\r
420//// with new I/O\r
421Int_t AliTPCclustererKr::FinderIO()\r
422{\r
423 // Krypton cluster finder for simulated events from MC\r
424\r
425 if (!fInput) { \r
426 Error("Digits2Clusters", "input tree not initialised");\r
427 return 10;\r
428 }\r
429 \r
430 if (!fOutput) {\r
431 Error("Digits2Clusters", "output tree not initialised");\r
432 return 11;\r
433 }\r
434\r
435 FindClusterKrIO();\r
436 return 0;\r
437}\r
438\r
98b76430 439\r
440\r
e2a4d72c 441Int_t AliTPCclustererKr::FinderIO(AliRawReader* rawReader)\r
442{\r
443 // Krypton cluster finder for the TPC raw data\r
98b76430 444 // this method is unsing AliAltroRawStreamV3\r
e2a4d72c 445 // fParam must be defined before\r
98b76430 446 \r
e2a4d72c 447 if(rawReader)fRawData=kTRUE; //set flag to data\r
98b76430 448 \r
e2a4d72c 449 if (!fOutput) {\r
450 Error("Digits2Clusters", "output tree not initialised");\r
451 return 11;\r
452 }\r
98b76430 453 \r
e2a4d72c 454 fParam->SetMaxTBin(fRecoParam->GetLastBin());//set number of timebins from reco -> param\r
455 // used later for memory allocation\r
456\r
98b76430 457 AliRawEventHeaderBase* eventHeader = (AliRawEventHeaderBase*)rawReader->GetEventHeader();\r
458 if (eventHeader){\r
459 fTimeStamp = eventHeader->Get("Timestamp");\r
f5178f60 460 fRun = rawReader->GetRunNumber();\r
98b76430 461 }\r
e2a4d72c 462\r
98b76430 463\r
464 Bool_t isAltro=kFALSE;\r
465 \r
e2a4d72c 466 AliTPCROC * roc = AliTPCROC::Instance();\r
467 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();\r
468 AliTPCAltroMapping** mapping =AliTPCcalibDB::Instance()->GetMapping();\r
469 //\r
98b76430 470 AliTPCRawStreamV3 input(rawReader,(AliAltroMapping**)mapping);\r
471 \r
e2a4d72c 472 const Int_t kNIS = fParam->GetNInnerSector();//number of inner sectors\r
473 const Int_t kNOS = fParam->GetNOuterSector();//number of outer sectors\r
474 const Int_t kNS = kNIS + kNOS;//all sectors\r
98b76430 475 \r
476 \r
477 //crate TPC view\r
478 AliTPCDigitsArray *digarr=new AliTPCDigitsArray(kFALSE);//data not sim\r
479 digarr->Setup(fParam);//as usually parameters\r
480 \r
481 for(Int_t iSec = 0; iSec < kNS; iSec++) {\r
482 AliTPCCalROC * noiseROC;\r
483 AliTPCCalROC noiseDummy(iSec);\r
484 if(noiseTPC==0x0){\r
485 noiseROC = &noiseDummy;//noise=0\r
486 }else{\r
487 noiseROC = noiseTPC->GetCalROC(iSec); // noise per given sector\r
488 }\r
489 Int_t nRows = 0; //number of rows in sector\r
490 Int_t nDDLs = 0; //number of DDLs\r
491 Int_t indexDDL = 0; //DDL index\r
492 if (iSec < kNIS) {\r
493 nRows = fParam->GetNRowLow();\r
494 nDDLs = 2;\r
495 indexDDL = iSec * 2;\r
496 }else {\r
497 nRows = fParam->GetNRowUp();\r
498 nDDLs = 4;\r
499 indexDDL = (iSec-kNIS) * 4 + kNIS * 2;\r
500 }\r
501 \r
502 //\r
503 // Load the raw data for corresponding DDLs\r
504 //\r
505 rawReader->Reset();\r
506 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);\r
507 \r
508 \r
509 while (input.NextDDL()){\r
510 // Allocate memory for rows in sector (pads(depends on row) x timebins)\r
511 if (!digarr->GetRow(iSec,0)){\r
512 for(Int_t iRow = 0; iRow < nRows; iRow++) {\r
513 digarr->CreateRow(iSec,iRow);\r
514 }//end loop over rows\r
515 }\r
516 //loop over pads\r
517 while ( input.NextChannel() ) {\r
518 Int_t iRow = input.GetRow();\r
519 Int_t iPad = input.GetPad();\r
520 //check row consistency\r
521 if (iRow < 0 ) continue;\r
522 if (iRow < 0 || iRow >= nRows){\r
523 AliError(Form("Pad-row index (%d) outside the range (%d -> %d) !",\r
524 iRow, 0, nRows -1));\r
525 continue;\r
526 }\r
527 \r
528 //check pad consistency\r
529 if (iPad < 0 || iPad >= (Int_t)(roc->GetNPads(iSec,iRow))) {\r
530 AliError(Form("Pad index (%d) outside the range (%d -> %d) !",\r
531 iPad, 0, roc->GetNPads(iSec,iRow) ));\r
532 continue;\r
533 }\r
534 \r
535 //loop over bunches\r
536 while ( input.NextBunch() ){\r
537 Int_t startTbin = (Int_t)input.GetStartTimeBin();\r
538 Int_t bunchlength = (Int_t)input.GetBunchLength();\r
539 const UShort_t *sig = input.GetSignals();\r
540 isAltro=kTRUE;\r
541 for (Int_t iTime = 0; iTime<bunchlength; iTime++){\r
542 Int_t iTimeBin=startTbin-iTime;\r
543 //\r
544 if(fDebugLevel==72){\r
545 fHistoRow->Fill(iRow);\r
546 fHistoPad->Fill(iPad);\r
547 fHistoTime->Fill(iTimeBin);\r
548 fHistoRowPad->Fill(iPad,iRow);\r
549 }else if(fDebugLevel>=0&&fDebugLevel<72){\r
550 if(iSec==fDebugLevel){\r
551 fHistoRow->Fill(iRow);\r
552 fHistoPad->Fill(iPad);\r
553 fHistoTime->Fill(iTimeBin);\r
554 fHistoRowPad->Fill(iPad,iRow);\r
555 }\r
556 }else if(fDebugLevel==73){\r
557 if(iSec<36){\r
558 fHistoRow->Fill(iRow);\r
559 fHistoPad->Fill(iPad);\r
560 fHistoTime->Fill(iTimeBin);\r
561 fHistoRowPad->Fill(iPad,iRow);\r
562 }\r
563 }else if(fDebugLevel==74){\r
564 if(iSec>=36){\r
565 fHistoRow->Fill(iRow);\r
566 fHistoPad->Fill(iPad);\r
567 fHistoTime->Fill(iTimeBin);\r
568 fHistoRowPad->Fill(iPad,iRow);\r
569 }\r
570 }\r
571 \r
572 //check time consistency\r
573 if ( iTimeBin < fRecoParam->GetFirstBin() || iTimeBin >= fRecoParam->GetLastBin()){\r
574 //cout<<iTimeBin<<endl;\r
575 continue;\r
576 AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",\r
577 iTimeBin, 0, fRecoParam->GetLastBin() -1));\r
578 }\r
579 //signal\r
580 Float_t signal=(Float_t)sig[iTime];\r
581 if (signal <= fZeroSup ||\r
582 iTimeBin < fFirstBin ||\r
583 iTimeBin > fLastBin\r
584 ) {\r
585 digarr->GetRow(iSec,iRow)->SetDigitFast(0,iTimeBin,iPad);\r
586 continue;\r
587 }\r
588 if (!noiseROC) continue;\r
589 Double_t noiseOnPad = noiseROC->GetValue(iRow,iPad);//noise on given pad and row in sector\r
590 if (noiseOnPad > fMaxNoiseAbs){\r
591 digarr->GetRow(iSec,iRow)->SetDigitFast(0,iTimeBin,iPad);\r
592 continue; // consider noisy pad as dead\r
593 }\r
594 if(signal <= fMaxNoiseSigma * noiseOnPad){\r
595 digarr->GetRow(iSec,iRow)->SetDigitFast(0,iTimeBin,iPad);\r
596 continue;\r
597 }\r
598 digarr->GetRow(iSec,iRow)->SetDigitFast(TMath::Nint(signal),iTimeBin,iPad);\r
599 }// end loop signals in bunch\r
600 }// end loop bunches\r
601 } // end loop pads\r
602 }// end ddl loop\r
603 }// end sector loop\r
604 SetDigArr(digarr);\r
605 if(isAltro) FindClusterKrIO();\r
606 delete digarr;\r
607 \r
608 return 0;\r
609}\r
e2a4d72c 610\r
611\r
98b76430 612\r
613\r
614\r
615Int_t AliTPCclustererKr::FinderIOold(AliRawReader* rawReader)\r
616{\r
617 // Krypton cluster finder for the TPC raw data\r
618 //\r
619 // fParam must be defined before\r
620 \r
621 if(rawReader)fRawData=kTRUE; //set flag to data\r
622 \r
623 if (!fOutput) {\r
624 Error("Digits2Clusters", "output tree not initialised");\r
625 return 11;\r
626 }\r
627 \r
628 fParam->SetMaxTBin(fRecoParam->GetLastBin());//set number of timebins from reco -> param\r
629 // used later for memory allocation\r
630\r
631 AliRawEventHeaderBase* eventHeader = (AliRawEventHeaderBase*)rawReader->GetEventHeader();\r
632 if (eventHeader){\r
633 fTimeStamp = eventHeader->Get("Timestamp");\r
f5178f60 634 fRun = rawReader->GetRunNumber();\r
98b76430 635 }\r
636 \r
637 Bool_t isAltro=kFALSE;\r
638 \r
639 AliTPCROC * roc = AliTPCROC::Instance();\r
640 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();\r
641 AliTPCAltroMapping** mapping =AliTPCcalibDB::Instance()->GetMapping();\r
642 //\r
643 AliTPCRawStream input(rawReader,(AliAltroMapping**)mapping);\r
644 \r
645 const Int_t kNIS = fParam->GetNInnerSector();//number of inner sectors\r
646 const Int_t kNOS = fParam->GetNOuterSector();//number of outer sectors\r
647 const Int_t kNS = kNIS + kNOS;//all sectors\r
648 \r
649 \r
e2a4d72c 650 //crate TPC view\r
651 AliTPCDigitsArray *digarr=new AliTPCDigitsArray(kFALSE);//data not sim\r
652 digarr->Setup(fParam);//as usually parameters\r
98b76430 653 \r
e2a4d72c 654 //\r
655 // Loop over sectors\r
656 //\r
657 for(Int_t iSec = 0; iSec < kNS; iSec++) {\r
658 AliTPCCalROC * noiseROC;\r
659 if(noiseTPC==0x0){\r
660 noiseROC =new AliTPCCalROC(iSec);//noise=0\r
661 }\r
662 else{\r
663 noiseROC = noiseTPC->GetCalROC(iSec); // noise per given sector\r
664 }\r
665 Int_t nRows = 0; //number of rows in sector\r
666 Int_t nDDLs = 0; //number of DDLs\r
667 Int_t indexDDL = 0; //DDL index\r
668 if (iSec < kNIS) {\r
669 nRows = fParam->GetNRowLow();\r
670 nDDLs = 2;\r
671 indexDDL = iSec * 2;\r
672 }else {\r
673 nRows = fParam->GetNRowUp();\r
674 nDDLs = 4;\r
675 indexDDL = (iSec-kNIS) * 4 + kNIS * 2;\r
676 }\r
98b76430 677 \r
e2a4d72c 678 //\r
679 // Load the raw data for corresponding DDLs\r
680 //\r
681 rawReader->Reset();\r
682 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);\r
98b76430 683 \r
e2a4d72c 684 if(input.Next()) {\r
685 isAltro=kTRUE;\r
686 // Allocate memory for rows in sector (pads(depends on row) x timebins)\r
687 for(Int_t iRow = 0; iRow < nRows; iRow++) {\r
98b76430 688 digarr->CreateRow(iSec,iRow);\r
e2a4d72c 689 }//end loop over rows\r
690 }\r
691 rawReader->Reset();\r
692 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);\r
98b76430 693 \r
e2a4d72c 694 //\r
695 // Begin loop over altro data\r
696 //\r
697 while (input.Next()) {\r
698 \r
699 //check sector consistency\r
700 if (input.GetSector() != iSec)\r
98b76430 701 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",iSec,input.GetSector()));\r
e2a4d72c 702 \r
703 Int_t iRow = input.GetRow();\r
704 Int_t iPad = input.GetPad();\r
705 Int_t iTimeBin = input.GetTime();\r
98b76430 706 \r
e2a4d72c 707 //\r
708 if(fDebugLevel==72){\r
98b76430 709 fHistoRow->Fill(iRow);\r
710 fHistoPad->Fill(iPad);\r
711 fHistoTime->Fill(iTimeBin);\r
712 fHistoRowPad->Fill(iPad,iRow);\r
e2a4d72c 713 }else if(fDebugLevel>=0&&fDebugLevel<72){\r
98b76430 714 if(iSec==fDebugLevel){\r
715 fHistoRow->Fill(iRow);\r
716 fHistoPad->Fill(iPad);\r
717 fHistoTime->Fill(iTimeBin);\r
718 fHistoRowPad->Fill(iPad,iRow);\r
719 }\r
e2a4d72c 720 }else if(fDebugLevel==73){\r
98b76430 721 if(iSec<36){\r
722 fHistoRow->Fill(iRow);\r
723 fHistoPad->Fill(iPad);\r
724 fHistoTime->Fill(iTimeBin);\r
725 fHistoRowPad->Fill(iPad,iRow);\r
726 }\r
e2a4d72c 727 }else if(fDebugLevel==74){\r
98b76430 728 if(iSec>=36){\r
729 fHistoRow->Fill(iRow);\r
730 fHistoPad->Fill(iPad);\r
731 fHistoTime->Fill(iTimeBin);\r
732 fHistoRowPad->Fill(iPad,iRow);\r
733 }\r
e2a4d72c 734 }\r
98b76430 735 \r
e2a4d72c 736 //check row consistency\r
737 if (iRow < 0 ) continue;\r
738 if (iRow < 0 || iRow >= nRows){\r
98b76430 739 AliError(Form("Pad-row index (%d) outside the range (%d -> %d) !",\r
740 iRow, 0, nRows -1));\r
741 continue;\r
e2a4d72c 742 }\r
98b76430 743 \r
e2a4d72c 744 //check pad consistency\r
745 if (iPad < 0 || iPad >= (Int_t)(roc->GetNPads(iSec,iRow))) {\r
98b76430 746 AliError(Form("Pad index (%d) outside the range (%d -> %d) !",\r
747 iPad, 0, roc->GetNPads(iSec,iRow) ));\r
748 continue;\r
e2a4d72c 749 }\r
98b76430 750 \r
e2a4d72c 751 //check time consistency\r
752 if ( iTimeBin < fRecoParam->GetFirstBin() || iTimeBin >= fRecoParam->GetLastBin()){\r
98b76430 753 //cout<<iTimeBin<<endl;\r
754 continue;\r
755 AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",\r
756 iTimeBin, 0, fRecoParam->GetLastBin() -1));\r
e2a4d72c 757 }\r
98b76430 758 \r
e2a4d72c 759 //signal\r
760 Int_t signal = input.GetSignal();\r
98b76430 761 if (signal <= fZeroSup ||\r
762 iTimeBin < fFirstBin ||\r
763 iTimeBin > fLastBin\r
764 ) {\r
765 digarr->GetRow(iSec,iRow)->SetDigitFast(0,iTimeBin,iPad);\r
766 continue;\r
767 }\r
e2a4d72c 768 if (!noiseROC) continue;\r
769 Double_t noiseOnPad = noiseROC->GetValue(iRow,iPad);//noise on given pad and row in sector\r
770 if (noiseOnPad > fMaxNoiseAbs){\r
98b76430 771 digarr->GetRow(iSec,iRow)->SetDigitFast(0,iTimeBin,iPad);\r
772 continue; // consider noisy pad as dead\r
e2a4d72c 773 }\r
774 if(signal <= fMaxNoiseSigma * noiseOnPad){\r
98b76430 775 digarr->GetRow(iSec,iRow)->SetDigitFast(0,iTimeBin,iPad);\r
776 continue;\r
e2a4d72c 777 }\r
98b76430 778 digarr->GetRow(iSec,iRow)->SetDigitFast(signal,iTimeBin,iPad);\r
e2a4d72c 779 }//end of loop over altro data\r
780 }//end of loop over sectors\r
781 \r
782 SetDigArr(digarr);\r
783 if(isAltro) FindClusterKrIO();\r
784 delete digarr;\r
98b76430 785 \r
e2a4d72c 786 return 0;\r
787}\r
788\r
789void AliTPCclustererKr::CleanSector(Int_t sector){\r
790 //\r
791 // clean isolated digits\r
792 // \r
793 const Int_t kNRows=fParam->GetNRow(sector);//number of rows in sector\r
794 for(Int_t iRow=0; iRow<kNRows; ++iRow){\r
795 AliSimDigits *digrow;\r
796 if(fRawData){\r
797 digrow = (AliSimDigits*)fDigarr->GetRow(sector,iRow);//real data\r
798 }else{\r
799 digrow = (AliSimDigits*)fDigarr->LoadRow(sector,iRow);//MC\r
800 }\r
801 if(!digrow) continue;\r
802 digrow->ExpandBuffer(); //decrunch\r
803 const Int_t kNPads = digrow->GetNCols(); // number of pads\r
804 const Int_t kNTime = digrow->GetNRows(); // number of timebins\r
805 for(Int_t iPad=1;iPad<kNPads-1;iPad++){\r
806 Short_t* val = digrow->GetDigitsColumn(iPad);\r
807\r
808 for(Int_t iTimeBin=1;iTimeBin<kNTime-1;iTimeBin++){\r
809 if (val[iTimeBin]<=0) continue;\r
810 if (val[iTimeBin-1]+val[iTimeBin+1]<fIsolCut) {val[iTimeBin]=0; continue;}\r
811 if (val[iTimeBin-kNTime]+val[iTimeBin+kNTime]<fIsolCut) {val[iTimeBin]=0; continue;}\r
812 //\r
813 if (val[iTimeBin-1-kNTime]+val[iTimeBin+1+kNTime]<fIsolCut) {val[iTimeBin]=0; continue;}\r
814 if (val[iTimeBin+1-kNTime]+val[iTimeBin-1+kNTime]<fIsolCut) {val[iTimeBin]=0; continue;}\r
815\r
816 }\r
817 }\r
818 }\r
819}\r
820\r
821\r
822////____________________________________________________________________________\r
823Int_t AliTPCclustererKr::FindClusterKrIO()\r
824{\r
825\r
826 //\r
827 //fParam and fDigarr must be set to run this method\r
828 //\r
829\r
830 Int_t clusterCounter=0;\r
831 const Int_t nTotalSector=fParam->GetNSector();//number of sectors\r
832 for(Int_t iSec=0; iSec<nTotalSector; ++iSec){\r
833 CleanSector(iSec);\r
834\r
835 //vector of maxima for each sector\r
836 //std::vector<AliPadMax*> maximaInSector;\r
837 TObjArray *maximaInSector=new TObjArray();//to store AliPadMax*\r
838\r
839 //\r
840 // looking for the maxima on the pad\r
841 //\r
842\r
843 const Int_t kNRows=fParam->GetNRow(iSec);//number of rows in sector\r
844 for(Int_t iRow=0; iRow<kNRows; ++iRow){\r
845 AliSimDigits *digrow;\r
846 if(fRawData){\r
847 digrow = (AliSimDigits*)fDigarr->GetRow(iSec,iRow);//real data\r
848 }else{\r
849 digrow = (AliSimDigits*)fDigarr->LoadRow(iSec,iRow);//MC\r
850 }\r
851 if(digrow){//if pointer exist\r
852 digrow->ExpandBuffer(); //decrunch\r
853 const Int_t kNPads = digrow->GetNCols(); // number of pads\r
854 const Int_t kNTime = digrow->GetNRows(); // number of timebins\r
855 for(Int_t iPad=0;iPad<kNPads;iPad++){\r
856 \r
857 Int_t timeBinMax=-1;//timebin of maximum \r
858 Int_t valueMaximum=-1;//value of maximum in adc\r
859 Int_t increaseBegin=-1;//timebin when increase starts\r
860 Int_t sumAdc=0;//sum of adc on the pad in maximum surrounding\r
861 bool ifIncreaseBegin=true;//flag - check if increasing started\r
862 bool ifMaximum=false;//flag - check if it could be maximum\r
863 Short_t* val = digrow->GetDigitsColumn(iPad);\r
864 for(Int_t iTimeBin=1;iTimeBin<kNTime-1;iTimeBin++){\r
865 if (!ifMaximum) {\r
866 if (val[iTimeBin]==-1) break; // 0 until the end\r
a6e0ebfe 867 for( ; iTimeBin<kNTime-2&&val[iTimeBin]<fMinAdc ;iTimeBin++) {}\r
e2a4d72c 868 }\r
869 //\r
870 Short_t adc = val[iTimeBin];\r
871\r
872 if(adc<fMinAdc){//standard was 3 for fMinAdc\r
873 if(ifMaximum){\r
874 if(iTimeBin-increaseBegin<fMinTimeBins){//at least 2 time bins\r
875 timeBinMax=-1;\r
876 valueMaximum=-1;\r
877 increaseBegin=-1;\r
878 sumAdc=0;\r
879 ifIncreaseBegin=true;\r
880 ifMaximum=false;\r
881 continue;\r
882 }\r
883 //insert maximum, default values and set flags\r
884 //Double_t xCord,yCord;\r
885 //GetXY(iSec,iRow,iPad,xCord,yCord);\r
886 Double_t x[]={iRow,iPad,iTimeBin};\r
887 Int_t i[]={iSec};\r
888 AliTPCTransform trafo;\r
889 trafo.Transform(x,i,0,1);\r
890 \r
891 AliPadMax *oneMaximum = new AliPadMax(AliTPCvtpr(valueMaximum,\r
892 timeBinMax,\r
893 iPad,\r
894 iRow,\r
895 x[0],//xCord,\r
896 x[1],//yCord,\r
897 x[2]/*timeBinMax*/),\r
898 increaseBegin,\r
899 iTimeBin-1,\r
900 sumAdc);\r
901 maximaInSector->AddLast(oneMaximum);\r
902 \r
903 timeBinMax=-1;\r
904 valueMaximum=-1;\r
905 increaseBegin=-1;\r
906 sumAdc=0;\r
907 ifIncreaseBegin=true;\r
908 ifMaximum=false;\r
909 }\r
910 continue;\r
911 }\r
912\r
913\r
914\r
915\r
916\r
917\r
918 if(ifIncreaseBegin){\r
919 ifIncreaseBegin=false;\r
920 increaseBegin=iTimeBin;\r
921 }\r
922 \r
923 if(adc>valueMaximum){\r
924 timeBinMax=iTimeBin;\r
925 valueMaximum=adc;\r
926 ifMaximum=true;\r
927 }\r
928 sumAdc+=adc;\r
929 if(iTimeBin==kNTime-1 && ifMaximum && kNTime-increaseBegin>fMinTimeBins){//on the edge\r
930 //at least 3 timebins\r
931 //insert maximum, default values and set flags\r
932 //Double_t xCord,yCord;\r
933 //GetXY(iSec,iRow,iPad,xCord,yCord);\r
934 Double_t x[]={iRow,iPad,iTimeBin};\r
935 Int_t i[]={iSec};\r
936 AliTPCTransform trafo;\r
937 trafo.Transform(x,i,0,1);\r
938 AliPadMax *oneMaximum = new AliPadMax(AliTPCvtpr(valueMaximum,\r
939 timeBinMax,\r
940 iPad,\r
941 iRow,\r
942 x[0],//xCord,\r
943 x[1],//yCord,\r
944 x[2]/*timeBinMax*/),\r
945 increaseBegin,\r
946 iTimeBin-1,\r
947 sumAdc);\r
948 maximaInSector->AddLast(oneMaximum);\r
949 \r
950 timeBinMax=-1;\r
951 valueMaximum=-1;\r
952 increaseBegin=-1;\r
953 sumAdc=0;\r
954 ifIncreaseBegin=true;\r
955 ifMaximum=false;\r
956 continue;\r
957 }\r
958 \r
959 }//end loop over timebins\r
960 }//end loop over pads\r
961// }else{\r
962// cout<<"Pointer does not exist!!"<<endl;\r
963 }//end if poiner exists\r
964 }//end loop over rows\r
965\r
966 MakeClusters(maximaInSector,iSec,clusterCounter);\r
967 //\r
968 maximaInSector->SetOwner(kTRUE);\r
969 maximaInSector->Delete();\r
970 delete maximaInSector;\r
971 }//end sector for\r
972 cout<<"Number of clusters in event: "<<clusterCounter<<endl;\r
973 return 0;\r
974}\r
975\r
976void AliTPCclustererKr::MakeClusters(TObjArray * maximaInSector, Int_t iSec, Int_t &clusterCounter){\r
977 //\r
978 // Make clusters\r
979 //\r
980\r
981 Int_t maxDig=0;\r
982 Int_t maxSumAdc=0;\r
983 Int_t maxTimeBin=0;\r
984 Int_t maxPad=0;\r
985 Int_t maxRow=0;\r
986 Double_t maxX=0;\r
987 Double_t maxY=0;\r
988 Double_t maxT=0;\r
989 Int_t entriesArr = maximaInSector->GetEntriesFast();\r
990 for(Int_t it1 = 0; it1 < entriesArr; ++it1 ) {\r
991 \r
992 AliPadMax *mp1=(AliPadMax *)maximaInSector->UncheckedAt(it1);\r
993 if (!mp1) continue;\r
994 AliTPCclusterKr clusterKr;\r
995 \r
996 Int_t nUsedPads=1;\r
997 Int_t clusterValue=0;\r
998 clusterValue+=(mp1)->GetSum();\r
999 list<Int_t> nUsedRows;\r
1000 nUsedRows.push_back((mp1)->GetRow());\r
1001 \r
1002 maxDig =(mp1)->GetAdc() ;\r
1003 maxSumAdc =(mp1)->GetSum() ;\r
1004 maxTimeBin =(mp1)->GetTime();\r
1005 maxPad =(mp1)->GetPad() ;\r
1006 maxRow =(mp1)->GetRow() ;\r
1007 maxX =(mp1)->GetX();\r
1008 maxY =(mp1)->GetY();\r
1009 maxT =(mp1)->GetT();\r
1010 \r
1011 AliSimDigits *digrowTmp;\r
1012 if(fRawData){\r
1013 digrowTmp = (AliSimDigits*)fDigarr->GetRow(iSec,(mp1)->GetRow());\r
1014 }else{\r
1015 digrowTmp = (AliSimDigits*)fDigarr->LoadRow(iSec,(mp1)->GetRow());\r
1016 }\r
1017 \r
1018 digrowTmp->ExpandBuffer(); //decrunch\r
1019 \r
1020 for(Int_t itb=(mp1)->GetBegin(); itb<((mp1)->GetEnd())+1; itb++){\r
1021 Int_t adcTmp = digrowTmp->GetDigitUnchecked(itb,(mp1)->GetPad());\r
1022 AliTPCvtpr *vtpr=new AliTPCvtpr(adcTmp,itb,(mp1)->GetPad(),(mp1)->GetRow(),(mp1)->GetX(),(mp1)->GetY(),(mp1)->GetT());\r
1023 clusterKr.AddDigitToCluster(vtpr);\r
1024 }\r
1025 clusterKr.SetCenter();//set centr of the cluster\r
1026 \r
1027 for(Int_t it2 = it1+1; it2 < entriesArr; ++it2 ) {\r
1028 AliPadMax *mp2=(AliPadMax *)maximaInSector->UncheckedAt(it2);\r
1029 if (!mp2) continue;\r
1030 if (TMath::Abs(clusterKr.GetCenterX() - (mp2)->GetX()) > fMaxPadRangeCm) continue;\r
1031 if (TMath::Abs(clusterKr.GetCenterY() - (mp2)->GetY()) > fMaxRowRangeCm) continue; \r
1032 if (TMath::Abs(clusterKr.GetCenterT() - (mp2)->GetT()) > fMaxTimeRange) continue;\r
1033\r
1034 {\r
1035 clusterValue+=(mp2)->GetSum();\r
1036 \r
1037 nUsedPads++;\r
1038 nUsedRows.push_back((mp2)->GetRow());\r
1039 \r
1040 AliSimDigits *digrowTmp1;\r
1041 if(fRawData){\r
1042 digrowTmp1 = (AliSimDigits*)fDigarr->GetRow(iSec,(mp2)->GetRow());\r
1043 }else{\r
1044 digrowTmp1 = (AliSimDigits*)fDigarr->LoadRow(iSec,(mp2)->GetRow());\r
1045 }\r
1046 digrowTmp1->ExpandBuffer(); //decrunch\r
1047 \r
1048 for(Int_t itb=(mp2)->GetBegin(); itb<(mp2)->GetEnd()+1; itb++){\r
1049 Int_t adcTmp = digrowTmp1->GetDigitUnchecked(itb,(mp2)->GetPad());\r
1050 AliTPCvtpr *vtpr=new AliTPCvtpr(adcTmp,itb,(mp2)->GetPad(),(mp2)->GetRow(),(mp2)->GetX(),(mp2)->GetY(),(mp2)->GetT());\r
1051 clusterKr.AddDigitToCluster(vtpr);\r
1052 }\r
1053 \r
1054 clusterKr.SetCenter();//set center of the cluster\r
1055 \r
1056 //which one is bigger\r
1057 if( (mp2)->GetAdc() > maxDig ){\r
1058 maxDig =(mp2)->GetAdc() ;\r
1059 maxSumAdc =(mp2)->GetSum() ;\r
1060 maxTimeBin =(mp2)->GetTime();\r
1061 maxPad =(mp2)->GetPad() ;\r
1062 maxRow =(mp2)->GetRow() ;\r
1063 maxX =(mp2)->GetX() ;\r
1064 maxY =(mp2)->GetY() ;\r
1065 maxT =(mp2)->GetT() ;\r
1066 } else if ( (mp2)->GetAdc() == maxDig ){\r
1067 if( (mp2)->GetSum() > maxSumAdc){\r
1068 maxDig =(mp2)->GetAdc() ;\r
1069 maxSumAdc =(mp2)->GetSum() ;\r
1070 maxTimeBin =(mp2)->GetTime();\r
1071 maxPad =(mp2)->GetPad() ;\r
1072 maxRow =(mp2)->GetRow() ;\r
1073 maxX =(mp2)->GetX() ;\r
1074 maxY =(mp2)->GetY() ;\r
1075 maxT =(mp2)->GetT() ;\r
1076 }\r
1077 }\r
1078 delete maximaInSector->RemoveAt(it2);\r
1079 }\r
1080 }//inside loop\r
1081 delete maximaInSector->RemoveAt(it1); \r
1082 clusterKr.SetSize();\r
1083 //through out clusters on the edge and noise\r
1084 //if(clusterValue/clusterKr.fCluster.size()<fValueToSize)continue;\r
1085 if(clusterValue/(clusterKr.GetSize())<fValueToSize)continue;\r
1086 \r
1087 clusterKr.SetADCcluster(clusterValue);\r
1088 clusterKr.SetNPads(nUsedPads);\r
1089 clusterKr.SetMax(AliTPCvtpr(maxDig,maxTimeBin,maxPad,maxRow,maxX,maxY,maxT));\r
1090 clusterKr.SetSec(iSec);\r
1091 clusterKr.SetSize();\r
1092 \r
1093 nUsedRows.sort();\r
1094 nUsedRows.unique();\r
1095 clusterKr.SetNRows(nUsedRows.size());\r
1096 clusterKr.SetCenter();\r
1097 \r
1098 clusterKr.SetRMS();//Set pad,row,timebin RMS\r
1099 clusterKr.Set1D();//Set size in pads and timebins\r
1100\r
98b76430 1101 clusterKr.SetTimeStamp(fTimeStamp);\r
f5178f60 1102 clusterKr.SetRun(fRun);\r
98b76430 1103\r
e2a4d72c 1104 clusterCounter++;\r
1105 \r
1106 \r
1107 //save each cluster into file\r
1108 if (fOutput){\r
1109 (*fOutput)<<"Kr"<<\r
1110 "Cl.="<<&clusterKr<<\r
1111 "\n";\r
1112 }\r
1113 //end of save each cluster into file adc.root\r
1114 }//outer loop\r
1115}\r
1116\r
1117\r
1118\r
1119////____________________________________________________________________________\r
1120\r
1121\r
1122void AliTPCclustererKr::GetXY(Int_t sec,Int_t row,Int_t pad,Double_t& xGlob,Double_t& yGlob){\r
1123 //\r
1124 //gives global XY coordinate of the pad\r
1125 //\r
1126\r
1127 Double_t yLocal = fParam->GetPadRowRadii(sec,row);//radius of row in sector in cm\r
1128\r
1129 Int_t padmax = fParam->GetNPads(sec,row);//number of pads in a given row\r
1130 Float_t padXSize;\r
1131 if(sec<fParam->GetNInnerSector())padXSize=0.4;\r
1132 else padXSize=0.6;\r
1133 Double_t xLocal=(pad+0.5-padmax/2.)*padXSize;//x-value of the center of pad\r
1134\r
1135 Float_t sin,cos;\r
1136 fParam->AdjustCosSin((Int_t)sec,cos,sin);//return sinus and cosinus of the sector\r
1137\r
1138 Double_t xGlob1 = xLocal * cos + yLocal * sin;\r
1139 Double_t yGlob1 = -xLocal * sin + yLocal * cos;\r
1140\r
1141\r
1142 Double_t rot=0;\r
1143 rot=TMath::Pi()/2.;\r
1144\r
1145 xGlob = xGlob1 * TMath::Cos(rot) + yGlob1 * TMath::Sin(rot);\r
1146 yGlob = -xGlob1 * TMath::Sin(rot) + yGlob1 * TMath::Cos(rot);\r
1147\r
1148 yGlob=-1*yGlob;\r
1149 if(sec<18||(sec>=36&&sec<54)) xGlob =-1*xGlob;\r
1150\r
1151\r
1152 return;\r
1153}\r