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