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