]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCclustererKr.cxx
Renamed
[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
831e26ce 620 if (!rawReader) return 1;\r
621\r
98b76430 622 if(rawReader)fRawData=kTRUE; //set flag to data\r
623 \r
624 if (!fOutput) {\r
625 Error("Digits2Clusters", "output tree not initialised");\r
626 return 11;\r
627 }\r
628 \r
629 fParam->SetMaxTBin(fRecoParam->GetLastBin());//set number of timebins from reco -> param\r
630 // used later for memory allocation\r
631\r
632 AliRawEventHeaderBase* eventHeader = (AliRawEventHeaderBase*)rawReader->GetEventHeader();\r
633 if (eventHeader){\r
634 fTimeStamp = eventHeader->Get("Timestamp");\r
f5178f60 635 fRun = rawReader->GetRunNumber();\r
98b76430 636 }\r
637 \r
638 Bool_t isAltro=kFALSE;\r
639 \r
640 AliTPCROC * roc = AliTPCROC::Instance();\r
641 AliTPCCalPad * noiseTPC = AliTPCcalibDB::Instance()->GetPadNoise();\r
642 AliTPCAltroMapping** mapping =AliTPCcalibDB::Instance()->GetMapping();\r
643 //\r
644 AliTPCRawStream input(rawReader,(AliAltroMapping**)mapping);\r
645 \r
646 const Int_t kNIS = fParam->GetNInnerSector();//number of inner sectors\r
647 const Int_t kNOS = fParam->GetNOuterSector();//number of outer sectors\r
648 const Int_t kNS = kNIS + kNOS;//all sectors\r
649 \r
650 \r
e2a4d72c 651 //crate TPC view\r
652 AliTPCDigitsArray *digarr=new AliTPCDigitsArray(kFALSE);//data not sim\r
653 digarr->Setup(fParam);//as usually parameters\r
98b76430 654 \r
e2a4d72c 655 //\r
656 // Loop over sectors\r
657 //\r
658 for(Int_t iSec = 0; iSec < kNS; iSec++) {\r
659 AliTPCCalROC * noiseROC;\r
660 if(noiseTPC==0x0){\r
661 noiseROC =new AliTPCCalROC(iSec);//noise=0\r
662 }\r
663 else{\r
664 noiseROC = noiseTPC->GetCalROC(iSec); // noise per given sector\r
665 }\r
666 Int_t nRows = 0; //number of rows in sector\r
667 Int_t nDDLs = 0; //number of DDLs\r
668 Int_t indexDDL = 0; //DDL index\r
669 if (iSec < kNIS) {\r
670 nRows = fParam->GetNRowLow();\r
671 nDDLs = 2;\r
672 indexDDL = iSec * 2;\r
673 }else {\r
674 nRows = fParam->GetNRowUp();\r
675 nDDLs = 4;\r
676 indexDDL = (iSec-kNIS) * 4 + kNIS * 2;\r
677 }\r
98b76430 678 \r
e2a4d72c 679 //\r
680 // Load the raw data for corresponding DDLs\r
681 //\r
682 rawReader->Reset();\r
683 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);\r
98b76430 684 \r
e2a4d72c 685 if(input.Next()) {\r
686 isAltro=kTRUE;\r
687 // Allocate memory for rows in sector (pads(depends on row) x timebins)\r
688 for(Int_t iRow = 0; iRow < nRows; iRow++) {\r
98b76430 689 digarr->CreateRow(iSec,iRow);\r
e2a4d72c 690 }//end loop over rows\r
691 }\r
692 rawReader->Reset();\r
693 rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);\r
98b76430 694 \r
e2a4d72c 695 //\r
696 // Begin loop over altro data\r
697 //\r
698 while (input.Next()) {\r
699 \r
700 //check sector consistency\r
701 if (input.GetSector() != iSec)\r
98b76430 702 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",iSec,input.GetSector()));\r
e2a4d72c 703 \r
704 Int_t iRow = input.GetRow();\r
705 Int_t iPad = input.GetPad();\r
706 Int_t iTimeBin = input.GetTime();\r
98b76430 707 \r
e2a4d72c 708 //\r
709 if(fDebugLevel==72){\r
98b76430 710 fHistoRow->Fill(iRow);\r
711 fHistoPad->Fill(iPad);\r
712 fHistoTime->Fill(iTimeBin);\r
713 fHistoRowPad->Fill(iPad,iRow);\r
e2a4d72c 714 }else if(fDebugLevel>=0&&fDebugLevel<72){\r
98b76430 715 if(iSec==fDebugLevel){\r
716 fHistoRow->Fill(iRow);\r
717 fHistoPad->Fill(iPad);\r
718 fHistoTime->Fill(iTimeBin);\r
719 fHistoRowPad->Fill(iPad,iRow);\r
720 }\r
e2a4d72c 721 }else if(fDebugLevel==73){\r
98b76430 722 if(iSec<36){\r
723 fHistoRow->Fill(iRow);\r
724 fHistoPad->Fill(iPad);\r
725 fHistoTime->Fill(iTimeBin);\r
726 fHistoRowPad->Fill(iPad,iRow);\r
727 }\r
e2a4d72c 728 }else if(fDebugLevel==74){\r
98b76430 729 if(iSec>=36){\r
730 fHistoRow->Fill(iRow);\r
731 fHistoPad->Fill(iPad);\r
732 fHistoTime->Fill(iTimeBin);\r
733 fHistoRowPad->Fill(iPad,iRow);\r
734 }\r
e2a4d72c 735 }\r
98b76430 736 \r
e2a4d72c 737 //check row consistency\r
738 if (iRow < 0 ) continue;\r
739 if (iRow < 0 || iRow >= nRows){\r
98b76430 740 AliError(Form("Pad-row index (%d) outside the range (%d -> %d) !",\r
741 iRow, 0, nRows -1));\r
742 continue;\r
e2a4d72c 743 }\r
98b76430 744 \r
e2a4d72c 745 //check pad consistency\r
746 if (iPad < 0 || iPad >= (Int_t)(roc->GetNPads(iSec,iRow))) {\r
98b76430 747 AliError(Form("Pad index (%d) outside the range (%d -> %d) !",\r
748 iPad, 0, roc->GetNPads(iSec,iRow) ));\r
749 continue;\r
e2a4d72c 750 }\r
98b76430 751 \r
e2a4d72c 752 //check time consistency\r
753 if ( iTimeBin < fRecoParam->GetFirstBin() || iTimeBin >= fRecoParam->GetLastBin()){\r
98b76430 754 //cout<<iTimeBin<<endl;\r
755 continue;\r
756 AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",\r
757 iTimeBin, 0, fRecoParam->GetLastBin() -1));\r
e2a4d72c 758 }\r
98b76430 759 \r
e2a4d72c 760 //signal\r
761 Int_t signal = input.GetSignal();\r
98b76430 762 if (signal <= fZeroSup ||\r
763 iTimeBin < fFirstBin ||\r
764 iTimeBin > fLastBin\r
765 ) {\r
766 digarr->GetRow(iSec,iRow)->SetDigitFast(0,iTimeBin,iPad);\r
767 continue;\r
768 }\r
e2a4d72c 769 if (!noiseROC) continue;\r
770 Double_t noiseOnPad = noiseROC->GetValue(iRow,iPad);//noise on given pad and row in sector\r
771 if (noiseOnPad > fMaxNoiseAbs){\r
98b76430 772 digarr->GetRow(iSec,iRow)->SetDigitFast(0,iTimeBin,iPad);\r
773 continue; // consider noisy pad as dead\r
e2a4d72c 774 }\r
775 if(signal <= fMaxNoiseSigma * noiseOnPad){\r
98b76430 776 digarr->GetRow(iSec,iRow)->SetDigitFast(0,iTimeBin,iPad);\r
777 continue;\r
e2a4d72c 778 }\r
98b76430 779 digarr->GetRow(iSec,iRow)->SetDigitFast(signal,iTimeBin,iPad);\r
e2a4d72c 780 }//end of loop over altro data\r
781 }//end of loop over sectors\r
782 \r
783 SetDigArr(digarr);\r
784 if(isAltro) FindClusterKrIO();\r
785 delete digarr;\r
98b76430 786 \r
e2a4d72c 787 return 0;\r
788}\r
789\r
790void AliTPCclustererKr::CleanSector(Int_t sector){\r
791 //\r
792 // clean isolated digits\r
793 // \r
794 const Int_t kNRows=fParam->GetNRow(sector);//number of rows in sector\r
795 for(Int_t iRow=0; iRow<kNRows; ++iRow){\r
796 AliSimDigits *digrow;\r
797 if(fRawData){\r
798 digrow = (AliSimDigits*)fDigarr->GetRow(sector,iRow);//real data\r
799 }else{\r
800 digrow = (AliSimDigits*)fDigarr->LoadRow(sector,iRow);//MC\r
801 }\r
802 if(!digrow) continue;\r
803 digrow->ExpandBuffer(); //decrunch\r
804 const Int_t kNPads = digrow->GetNCols(); // number of pads\r
805 const Int_t kNTime = digrow->GetNRows(); // number of timebins\r
806 for(Int_t iPad=1;iPad<kNPads-1;iPad++){\r
807 Short_t* val = digrow->GetDigitsColumn(iPad);\r
808\r
809 for(Int_t iTimeBin=1;iTimeBin<kNTime-1;iTimeBin++){\r
810 if (val[iTimeBin]<=0) continue;\r
811 if (val[iTimeBin-1]+val[iTimeBin+1]<fIsolCut) {val[iTimeBin]=0; continue;}\r
812 if (val[iTimeBin-kNTime]+val[iTimeBin+kNTime]<fIsolCut) {val[iTimeBin]=0; continue;}\r
813 //\r
814 if (val[iTimeBin-1-kNTime]+val[iTimeBin+1+kNTime]<fIsolCut) {val[iTimeBin]=0; continue;}\r
815 if (val[iTimeBin+1-kNTime]+val[iTimeBin-1+kNTime]<fIsolCut) {val[iTimeBin]=0; continue;}\r
816\r
817 }\r
818 }\r
819 }\r
820}\r
821\r
822\r
823////____________________________________________________________________________\r
824Int_t AliTPCclustererKr::FindClusterKrIO()\r
825{\r
826\r
827 //\r
828 //fParam and fDigarr must be set to run this method\r
829 //\r
830\r
831 Int_t clusterCounter=0;\r
832 const Int_t nTotalSector=fParam->GetNSector();//number of sectors\r
833 for(Int_t iSec=0; iSec<nTotalSector; ++iSec){\r
834 CleanSector(iSec);\r
835\r
836 //vector of maxima for each sector\r
837 //std::vector<AliPadMax*> maximaInSector;\r
838 TObjArray *maximaInSector=new TObjArray();//to store AliPadMax*\r
839\r
840 //\r
841 // looking for the maxima on the pad\r
842 //\r
843\r
844 const Int_t kNRows=fParam->GetNRow(iSec);//number of rows in sector\r
845 for(Int_t iRow=0; iRow<kNRows; ++iRow){\r
846 AliSimDigits *digrow;\r
847 if(fRawData){\r
848 digrow = (AliSimDigits*)fDigarr->GetRow(iSec,iRow);//real data\r
849 }else{\r
850 digrow = (AliSimDigits*)fDigarr->LoadRow(iSec,iRow);//MC\r
851 }\r
852 if(digrow){//if pointer exist\r
853 digrow->ExpandBuffer(); //decrunch\r
854 const Int_t kNPads = digrow->GetNCols(); // number of pads\r
855 const Int_t kNTime = digrow->GetNRows(); // number of timebins\r
856 for(Int_t iPad=0;iPad<kNPads;iPad++){\r
857 \r
858 Int_t timeBinMax=-1;//timebin of maximum \r
859 Int_t valueMaximum=-1;//value of maximum in adc\r
860 Int_t increaseBegin=-1;//timebin when increase starts\r
861 Int_t sumAdc=0;//sum of adc on the pad in maximum surrounding\r
862 bool ifIncreaseBegin=true;//flag - check if increasing started\r
863 bool ifMaximum=false;//flag - check if it could be maximum\r
864 Short_t* val = digrow->GetDigitsColumn(iPad);\r
865 for(Int_t iTimeBin=1;iTimeBin<kNTime-1;iTimeBin++){\r
866 if (!ifMaximum) {\r
867 if (val[iTimeBin]==-1) break; // 0 until the end\r
a6e0ebfe 868 for( ; iTimeBin<kNTime-2&&val[iTimeBin]<fMinAdc ;iTimeBin++) {}\r
e2a4d72c 869 }\r
870 //\r
871 Short_t adc = val[iTimeBin];\r
872\r
873 if(adc<fMinAdc){//standard was 3 for fMinAdc\r
874 if(ifMaximum){\r
875 if(iTimeBin-increaseBegin<fMinTimeBins){//at least 2 time bins\r
876 timeBinMax=-1;\r
877 valueMaximum=-1;\r
878 increaseBegin=-1;\r
879 sumAdc=0;\r
880 ifIncreaseBegin=true;\r
881 ifMaximum=false;\r
882 continue;\r
883 }\r
884 //insert maximum, default values and set flags\r
885 //Double_t xCord,yCord;\r
886 //GetXY(iSec,iRow,iPad,xCord,yCord);\r
887 Double_t x[]={iRow,iPad,iTimeBin};\r
888 Int_t i[]={iSec};\r
435a15ae 889 AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;\r
890\r
891 transform->Transform(x,i,0,1);\r
e2a4d72c 892 \r
893 AliPadMax *oneMaximum = new AliPadMax(AliTPCvtpr(valueMaximum,\r
894 timeBinMax,\r
895 iPad,\r
896 iRow,\r
897 x[0],//xCord,\r
898 x[1],//yCord,\r
899 x[2]/*timeBinMax*/),\r
900 increaseBegin,\r
901 iTimeBin-1,\r
902 sumAdc);\r
903 maximaInSector->AddLast(oneMaximum);\r
904 \r
905 timeBinMax=-1;\r
906 valueMaximum=-1;\r
907 increaseBegin=-1;\r
908 sumAdc=0;\r
909 ifIncreaseBegin=true;\r
910 ifMaximum=false;\r
911 }\r
912 continue;\r
913 }\r
914\r
915\r
916\r
917\r
918\r
919\r
920 if(ifIncreaseBegin){\r
921 ifIncreaseBegin=false;\r
922 increaseBegin=iTimeBin;\r
923 }\r
924 \r
925 if(adc>valueMaximum){\r
926 timeBinMax=iTimeBin;\r
927 valueMaximum=adc;\r
928 ifMaximum=true;\r
929 }\r
930 sumAdc+=adc;\r
931 if(iTimeBin==kNTime-1 && ifMaximum && kNTime-increaseBegin>fMinTimeBins){//on the edge\r
932 //at least 3 timebins\r
933 //insert maximum, default values and set flags\r
934 //Double_t xCord,yCord;\r
935 //GetXY(iSec,iRow,iPad,xCord,yCord);\r
936 Double_t x[]={iRow,iPad,iTimeBin};\r
937 Int_t i[]={iSec};\r
435a15ae 938 //AliTPCTransform trafo;\r
939 //trafo.Transform(x,i,0,1);\r
940\r
941 AliTPCTransform *transform = AliTPCcalibDB::Instance()->GetTransform() ;\r
942\r
943 transform->Transform(x,i,0,1);\r
944\r
e2a4d72c 945 AliPadMax *oneMaximum = new AliPadMax(AliTPCvtpr(valueMaximum,\r
946 timeBinMax,\r
947 iPad,\r
948 iRow,\r
949 x[0],//xCord,\r
950 x[1],//yCord,\r
951 x[2]/*timeBinMax*/),\r
952 increaseBegin,\r
953 iTimeBin-1,\r
954 sumAdc);\r
955 maximaInSector->AddLast(oneMaximum);\r
956 \r
957 timeBinMax=-1;\r
958 valueMaximum=-1;\r
959 increaseBegin=-1;\r
960 sumAdc=0;\r
961 ifIncreaseBegin=true;\r
962 ifMaximum=false;\r
963 continue;\r
964 }\r
965 \r
966 }//end loop over timebins\r
967 }//end loop over pads\r
968// }else{\r
969// cout<<"Pointer does not exist!!"<<endl;\r
970 }//end if poiner exists\r
971 }//end loop over rows\r
972\r
973 MakeClusters(maximaInSector,iSec,clusterCounter);\r
974 //\r
975 maximaInSector->SetOwner(kTRUE);\r
976 maximaInSector->Delete();\r
977 delete maximaInSector;\r
978 }//end sector for\r
979 cout<<"Number of clusters in event: "<<clusterCounter<<endl;\r
980 return 0;\r
981}\r
982\r
983void AliTPCclustererKr::MakeClusters(TObjArray * maximaInSector, Int_t iSec, Int_t &clusterCounter){\r
984 //\r
985 // Make clusters\r
986 //\r
987\r
988 Int_t maxDig=0;\r
989 Int_t maxSumAdc=0;\r
990 Int_t maxTimeBin=0;\r
991 Int_t maxPad=0;\r
992 Int_t maxRow=0;\r
993 Double_t maxX=0;\r
994 Double_t maxY=0;\r
995 Double_t maxT=0;\r
996 Int_t entriesArr = maximaInSector->GetEntriesFast();\r
997 for(Int_t it1 = 0; it1 < entriesArr; ++it1 ) {\r
998 \r
999 AliPadMax *mp1=(AliPadMax *)maximaInSector->UncheckedAt(it1);\r
1000 if (!mp1) continue;\r
1001 AliTPCclusterKr clusterKr;\r
1002 \r
1003 Int_t nUsedPads=1;\r
1004 Int_t clusterValue=0;\r
1005 clusterValue+=(mp1)->GetSum();\r
1006 list<Int_t> nUsedRows;\r
1007 nUsedRows.push_back((mp1)->GetRow());\r
1008 \r
1009 maxDig =(mp1)->GetAdc() ;\r
1010 maxSumAdc =(mp1)->GetSum() ;\r
1011 maxTimeBin =(mp1)->GetTime();\r
1012 maxPad =(mp1)->GetPad() ;\r
1013 maxRow =(mp1)->GetRow() ;\r
1014 maxX =(mp1)->GetX();\r
1015 maxY =(mp1)->GetY();\r
1016 maxT =(mp1)->GetT();\r
1017 \r
1018 AliSimDigits *digrowTmp;\r
1019 if(fRawData){\r
1020 digrowTmp = (AliSimDigits*)fDigarr->GetRow(iSec,(mp1)->GetRow());\r
1021 }else{\r
1022 digrowTmp = (AliSimDigits*)fDigarr->LoadRow(iSec,(mp1)->GetRow());\r
1023 }\r
1024 \r
1025 digrowTmp->ExpandBuffer(); //decrunch\r
1026 \r
1027 for(Int_t itb=(mp1)->GetBegin(); itb<((mp1)->GetEnd())+1; itb++){\r
1028 Int_t adcTmp = digrowTmp->GetDigitUnchecked(itb,(mp1)->GetPad());\r
1029 AliTPCvtpr *vtpr=new AliTPCvtpr(adcTmp,itb,(mp1)->GetPad(),(mp1)->GetRow(),(mp1)->GetX(),(mp1)->GetY(),(mp1)->GetT());\r
1030 clusterKr.AddDigitToCluster(vtpr);\r
1031 }\r
1032 clusterKr.SetCenter();//set centr of the cluster\r
1033 \r
1034 for(Int_t it2 = it1+1; it2 < entriesArr; ++it2 ) {\r
1035 AliPadMax *mp2=(AliPadMax *)maximaInSector->UncheckedAt(it2);\r
1036 if (!mp2) continue;\r
1037 if (TMath::Abs(clusterKr.GetCenterX() - (mp2)->GetX()) > fMaxPadRangeCm) continue;\r
1038 if (TMath::Abs(clusterKr.GetCenterY() - (mp2)->GetY()) > fMaxRowRangeCm) continue; \r
1039 if (TMath::Abs(clusterKr.GetCenterT() - (mp2)->GetT()) > fMaxTimeRange) continue;\r
1040\r
1041 {\r
1042 clusterValue+=(mp2)->GetSum();\r
1043 \r
1044 nUsedPads++;\r
1045 nUsedRows.push_back((mp2)->GetRow());\r
1046 \r
1047 AliSimDigits *digrowTmp1;\r
1048 if(fRawData){\r
1049 digrowTmp1 = (AliSimDigits*)fDigarr->GetRow(iSec,(mp2)->GetRow());\r
1050 }else{\r
1051 digrowTmp1 = (AliSimDigits*)fDigarr->LoadRow(iSec,(mp2)->GetRow());\r
1052 }\r
1053 digrowTmp1->ExpandBuffer(); //decrunch\r
1054 \r
1055 for(Int_t itb=(mp2)->GetBegin(); itb<(mp2)->GetEnd()+1; itb++){\r
1056 Int_t adcTmp = digrowTmp1->GetDigitUnchecked(itb,(mp2)->GetPad());\r
1057 AliTPCvtpr *vtpr=new AliTPCvtpr(adcTmp,itb,(mp2)->GetPad(),(mp2)->GetRow(),(mp2)->GetX(),(mp2)->GetY(),(mp2)->GetT());\r
1058 clusterKr.AddDigitToCluster(vtpr);\r
1059 }\r
1060 \r
1061 clusterKr.SetCenter();//set center of the cluster\r
1062 \r
1063 //which one is bigger\r
1064 if( (mp2)->GetAdc() > maxDig ){\r
1065 maxDig =(mp2)->GetAdc() ;\r
1066 maxSumAdc =(mp2)->GetSum() ;\r
1067 maxTimeBin =(mp2)->GetTime();\r
1068 maxPad =(mp2)->GetPad() ;\r
1069 maxRow =(mp2)->GetRow() ;\r
1070 maxX =(mp2)->GetX() ;\r
1071 maxY =(mp2)->GetY() ;\r
1072 maxT =(mp2)->GetT() ;\r
1073 } else if ( (mp2)->GetAdc() == maxDig ){\r
1074 if( (mp2)->GetSum() > maxSumAdc){\r
1075 maxDig =(mp2)->GetAdc() ;\r
1076 maxSumAdc =(mp2)->GetSum() ;\r
1077 maxTimeBin =(mp2)->GetTime();\r
1078 maxPad =(mp2)->GetPad() ;\r
1079 maxRow =(mp2)->GetRow() ;\r
1080 maxX =(mp2)->GetX() ;\r
1081 maxY =(mp2)->GetY() ;\r
1082 maxT =(mp2)->GetT() ;\r
1083 }\r
1084 }\r
1085 delete maximaInSector->RemoveAt(it2);\r
1086 }\r
1087 }//inside loop\r
1088 delete maximaInSector->RemoveAt(it1); \r
1089 clusterKr.SetSize();\r
1090 //through out clusters on the edge and noise\r
1091 //if(clusterValue/clusterKr.fCluster.size()<fValueToSize)continue;\r
1092 if(clusterValue/(clusterKr.GetSize())<fValueToSize)continue;\r
1093 \r
1094 clusterKr.SetADCcluster(clusterValue);\r
1095 clusterKr.SetNPads(nUsedPads);\r
1096 clusterKr.SetMax(AliTPCvtpr(maxDig,maxTimeBin,maxPad,maxRow,maxX,maxY,maxT));\r
1097 clusterKr.SetSec(iSec);\r
1098 clusterKr.SetSize();\r
1099 \r
1100 nUsedRows.sort();\r
1101 nUsedRows.unique();\r
1102 clusterKr.SetNRows(nUsedRows.size());\r
1103 clusterKr.SetCenter();\r
1104 \r
1105 clusterKr.SetRMS();//Set pad,row,timebin RMS\r
1106 clusterKr.Set1D();//Set size in pads and timebins\r
1107\r
98b76430 1108 clusterKr.SetTimeStamp(fTimeStamp);\r
f5178f60 1109 clusterKr.SetRun(fRun);\r
98b76430 1110\r
e2a4d72c 1111 clusterCounter++;\r
1112 \r
1113 \r
1114 //save each cluster into file\r
1115 if (fOutput){\r
1116 (*fOutput)<<"Kr"<<\r
1117 "Cl.="<<&clusterKr<<\r
1118 "\n";\r
1119 }\r
1120 //end of save each cluster into file adc.root\r
1121 }//outer loop\r
1122}\r
1123\r
1124\r
1125\r
1126////____________________________________________________________________________\r
1127\r
1128\r
1129void AliTPCclustererKr::GetXY(Int_t sec,Int_t row,Int_t pad,Double_t& xGlob,Double_t& yGlob){\r
1130 //\r
1131 //gives global XY coordinate of the pad\r
1132 //\r
1133\r
1134 Double_t yLocal = fParam->GetPadRowRadii(sec,row);//radius of row in sector in cm\r
1135\r
1136 Int_t padmax = fParam->GetNPads(sec,row);//number of pads in a given row\r
1137 Float_t padXSize;\r
1138 if(sec<fParam->GetNInnerSector())padXSize=0.4;\r
1139 else padXSize=0.6;\r
1140 Double_t xLocal=(pad+0.5-padmax/2.)*padXSize;//x-value of the center of pad\r
1141\r
1142 Float_t sin,cos;\r
1143 fParam->AdjustCosSin((Int_t)sec,cos,sin);//return sinus and cosinus of the sector\r
1144\r
1145 Double_t xGlob1 = xLocal * cos + yLocal * sin;\r
1146 Double_t yGlob1 = -xLocal * sin + yLocal * cos;\r
1147\r
1148\r
1149 Double_t rot=0;\r
1150 rot=TMath::Pi()/2.;\r
1151\r
1152 xGlob = xGlob1 * TMath::Cos(rot) + yGlob1 * TMath::Sin(rot);\r
1153 yGlob = -xGlob1 * TMath::Sin(rot) + yGlob1 * TMath::Cos(rot);\r
1154\r
1155 yGlob=-1*yGlob;\r
1156 if(sec<18||(sec>=36&&sec<54)) xGlob =-1*xGlob;\r
1157\r
1158\r
1159 return;\r
1160}\r