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