]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PMD/AliPMDClusterFinder.cxx
correcting baryon mass subtraction for visible energy in MC
[u/mrichter/AliRoot.git] / PMD / AliPMDClusterFinder.cxx
CommitLineData
ed228cbc 1/***************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
15
01709453 16//-----------------------------------------------------//
17// //
18// Date : August 05 2003 //
19// This reads the file PMD.digits.root(TreeD), //
20// calls the Clustering algorithm and stores the //
21// clustering output in PMD.RecPoints.root(TreeR) //
22// //
23//-----------------------------------------------------//
24
25#include <Riostream.h>
01709453 26#include <TTree.h>
01709453 27#include <TObjArray.h>
28#include <TClonesArray.h>
08cf9354 29#include <TSystem.h>
01709453 30
7e9508a7 31#include "AliLog.h"
01709453 32#include "AliRunLoader.h"
33#include "AliLoader.h"
5f55af10 34#include "AliRawReader.h"
01709453 35
36#include "AliPMDdigit.h"
37#include "AliPMDClusterFinder.h"
38#include "AliPMDClustering.h"
5c5cadd5 39#include "AliPMDClusteringV1.h"
01709453 40#include "AliPMDcluster.h"
96377d57 41#include "AliPMDrecpoint1.h"
5c5cadd5 42#include "AliPMDrechit.h"
5f55af10 43#include "AliPMDRawStream.h"
09a06455 44#include "AliPMDCalibData.h"
35535af7 45#include "AliPMDPedestal.h"
df42ab21 46#include "AliPMDddldata.h"
7adbd8e9 47#include "AliPMDHotData.h"
963fe679 48#include "AliPMDNoiseCut.h"
b9e2f33c 49#include "AliPMDRecoParam.h"
50#include "AliPMDReconstructor.h"
09a06455 51
2332574a 52#include "AliDAQ.h"
09a06455 53#include "AliCDBManager.h"
54#include "AliCDBEntry.h"
01709453 55
7e9508a7 56
57
01709453 58ClassImp(AliPMDClusterFinder)
b208c6a3 59
7e9508a7 60AliPMDClusterFinder::AliPMDClusterFinder():
61 fRunLoader(0),
62 fPMDLoader(0),
35535af7 63 fCalibGain(GetCalibGain()),
64 fCalibPed(GetCalibPed()),
7adbd8e9 65 fCalibHot(GetCalibHot()),
963fe679 66 fNoiseCut(GetNoiseCut()),
b9e2f33c 67 fRecoParam(0x0),
7e9508a7 68 fTreeD(0),
69 fTreeR(0),
70 fDigits(new TClonesArray("AliPMDdigit", 1000)),
71 fRecpoints(new TClonesArray("AliPMDrecpoint1", 1000)),
5c5cadd5 72 fRechits(new TClonesArray("AliPMDrechit", 1000)),
7e9508a7 73 fNpoint(0),
5c5cadd5 74 fNhit(0),
963fe679 75 fDetNo(0)
7e9508a7 76{
77//
78// Constructor
79//
80}
81// ------------------------------------------------------------------------- //
dfaeee5f 82AliPMDClusterFinder::AliPMDClusterFinder(AliRunLoader* runLoader):
83 fRunLoader(runLoader),
84 fPMDLoader(runLoader->GetLoader("PMDLoader")),
35535af7 85 fCalibGain(GetCalibGain()),
86 fCalibPed(GetCalibPed()),
7adbd8e9 87 fCalibHot(GetCalibHot()),
963fe679 88 fNoiseCut(GetNoiseCut()),
b9e2f33c 89 fRecoParam(0x0),
1758e4fe 90 fTreeD(0),
91 fTreeR(0),
ebd83c56 92 fDigits(new TClonesArray("AliPMDdigit", 1000)),
1758e4fe 93 fRecpoints(new TClonesArray("AliPMDrecpoint1", 1000)),
5c5cadd5 94 fRechits(new TClonesArray("AliPMDrechit", 1000)),
1758e4fe 95 fNpoint(0),
5c5cadd5 96 fNhit(0),
963fe679 97 fDetNo(0)
01709453 98{
b208c6a3 99//
dfaeee5f 100// Constructor
b208c6a3 101//
a48edddd 102}
103// ------------------------------------------------------------------------- //
2332574a 104AliPMDClusterFinder::AliPMDClusterFinder(const AliPMDClusterFinder & finder):
105 TObject(finder),
106 fRunLoader(0),
107 fPMDLoader(0),
35535af7 108 fCalibGain(GetCalibGain()),
109 fCalibPed(GetCalibPed()),
7adbd8e9 110 fCalibHot(GetCalibHot()),
963fe679 111 fNoiseCut(GetNoiseCut()),
b9e2f33c 112 fRecoParam(0x0),
2332574a 113 fTreeD(0),
114 fTreeR(0),
115 fDigits(NULL),
116 fRecpoints(NULL),
117 fRechits(NULL),
118 fNpoint(0),
119 fNhit(0),
963fe679 120 fDetNo(0)
a48edddd 121{
122 // copy constructor
123 AliError("Copy constructor not allowed");
124}
125// ------------------------------------------------------------------------- //
126AliPMDClusterFinder &AliPMDClusterFinder::operator=(const AliPMDClusterFinder & /*finder*/)
127{
128 // assignment op
129 AliError("Assignment Operator not allowed");
130 return *this;
01709453 131}
1758e4fe 132// ------------------------------------------------------------------------- //
01709453 133AliPMDClusterFinder::~AliPMDClusterFinder()
134{
b208c6a3 135 // Destructor
ebd83c56 136 if (fDigits)
137 {
b9e2f33c 138 fDigits->Clear();
ebd83c56 139 }
1758e4fe 140 if (fRecpoints)
141 {
722ccc67 142 fRecpoints->Clear();
1758e4fe 143 }
5c5cadd5 144 if (fRechits)
145 {
722ccc67 146 fRechits->Clear();
5c5cadd5 147 }
5254fcca 148
01709453 149}
1758e4fe 150// ------------------------------------------------------------------------- //
b208c6a3 151
8fbad6d3 152void AliPMDClusterFinder::Digits2RecPoints(TTree *digitsTree,
153 TTree *clustersTree)
154{
155 // Converts digits to recpoints after running clustering
156 // algorithm on CPV plane and PREshower plane
2c1131dd 157 //
f0290c15 158 // This algorithm is called during the reconstruction from digits
8fbad6d3 159
be8b7039 160 Int_t det = 0, smn = 0;
161 Int_t xpos = 0, ypos = 0;
162 Int_t ismn = 0;
163 Int_t idet = 0;
164 Float_t adc = 0.;
165 Float_t clusdata[6] = {0.,0.,0.,0.,0.,0.};
8fbad6d3 166
5254fcca 167 AliPMDcluster *pmdcl = 0x0;
168
8fbad6d3 169 TObjArray *pmdcont = new TObjArray();
27a78deb 170
8fbad6d3 171 AliPMDClustering *pmdclust = new AliPMDClusteringV1();
172
b9e2f33c 173 // Fetch the reco param object
174
175 fRecoParam = AliPMDReconstructor::GetRecoParam();
176 if(fRecoParam == 0x0)
177 {
178 AliFatal("No Reco Param found for PMD!!!");
179 }
180
8fbad6d3 181
8fbad6d3 182 AliPMDdigit *pmddigit;
183 TBranch *branch = digitsTree->GetBranch("PMDDigit");
184 branch->SetAddress(&fDigits);
185
186 ResetRecpoint();
187
8fbad6d3 188 Int_t bufsize = 16000;
189 TBranch * branch1 = clustersTree->Branch("PMDRecpoint", &fRecpoints, bufsize);
190 TBranch * branch2 = clustersTree->Branch("PMDRechit", &fRechits, bufsize);
191
192 Int_t nmodules = (Int_t) digitsTree->GetEntries();
193
194 for (Int_t imodule = 0; imodule < nmodules; imodule++)
195 {
5254fcca 196
197 Int_t totADCMod = 0;
8fbad6d3 198 ResetCellADC();
199 digitsTree->GetEntry(imodule);
200 Int_t nentries = fDigits->GetLast();
201 for (Int_t ient = 0; ient < nentries+1; ient++)
202 {
203 pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
204
205 det = pmddigit->GetDetector();
206 smn = pmddigit->GetSMNumber();
207 xpos = pmddigit->GetRow();
208 ypos = pmddigit->GetColumn();
209 adc = pmddigit->GetADC();
8a523d22 210
211 if(det < 0 || det > 1)
212 {
213 AliError(Form("*CPV/PRE NUMBER WRONG %d *",det));
214 continue;
215 }
72f10132 216 if(smn == -1 || smn > 23)
8a523d22 217 {
218 AliError(Form("*MODULE NUMBER WRONG %d *",smn));
219 continue;
220 }
221
72f10132 222 if(xpos < 0 || xpos > 47 || ypos < 0 || ypos > 95)
f0290c15 223 {
224 AliError(Form("*Row %d and Column NUMBER %d NOT Valid *",
225 xpos, ypos));
226 continue;
227 }
228
35535af7 229 // Pedestal Subtraction
230 Int_t pedmeanrms = fCalibPed->GetPedMeanRms(det,smn,xpos,ypos);
37edc588 231 Int_t pedrms1 = (Int_t) pedmeanrms%100;
35535af7 232 Float_t pedrms = (Float_t)pedrms1/10.;
233 Float_t pedmean = (Float_t) (pedmeanrms - pedrms1)/1000.0;
234 //printf("%f %f\n",pedmean, pedrms);
235
236 Float_t adc1 = adc - (pedmean + 3.0*pedrms);
237
7adbd8e9 238 // Hot cell - set the cell adc = 0
239 Float_t hotflag = fCalibHot->GetHotChannel(det,smn,xpos,ypos);
8a523d22 240 if (hotflag == 1.) adc1 = 0;
7adbd8e9 241
8fbad6d3 242 // CALIBRATION
35535af7 243 Float_t gain = fCalibGain->GetGainFact(det,smn,xpos,ypos);
8fbad6d3 244 // printf("adc = %d gain = %f\n",adc,gain);
35535af7 245
246 adc = adc1*gain;
8fbad6d3 247
920e13db 248 fCellTrack[xpos][ypos] = pmddigit->GetTrackNumber();
249 fCellPid[xpos][ypos] = pmddigit->GetTrackPid();
8fbad6d3 250 fCellADC[xpos][ypos] = (Double_t) adc;
5254fcca 251
5fbd2db1 252 totADCMod += (Int_t) adc;
5254fcca 253
8fbad6d3 254 }
255
256 idet = det;
257 ismn = smn;
5254fcca 258
259 if (totADCMod <= 0) continue;
260
b9e2f33c 261 // Set the minimum noise cut per module before clustering
262
963fe679 263 // Int_t imod = idet*24 + ismn;
264
265
266 // Int_t cluspar = fRecoParam->GetPbPbParam()->GetClusteringParam();
4755c3f1 267 Int_t cluspar = fRecoParam->GetPPParam()->GetClusteringParam();
963fe679 268 // Int_t cluspar = fRecoParam->GetCosmicParam()->GetClusteringParam();
b9e2f33c 269
4755c3f1 270 pmdclust->SetClusteringParam(cluspar);
b9e2f33c 271
4e642bdd 272 Float_t encut = 4.;
963fe679 273 pmdclust->SetEdepCut(encut);
27a78deb 274 pmdclust->DoClust(idet,ismn,fCellTrack,fCellPid,fCellADC,pmdcont);
8fbad6d3 275
276 Int_t nentries1 = pmdcont->GetEntries();
277
278 AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
279
280 for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
281 {
920e13db 282 pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
8fbad6d3 283 idet = pmdcl->GetDetector();
284 ismn = pmdcl->GetSMN();
285 clusdata[0] = pmdcl->GetClusX();
286 clusdata[1] = pmdcl->GetClusY();
287 clusdata[2] = pmdcl->GetClusADC();
288 clusdata[3] = pmdcl->GetClusCells();
289 clusdata[4] = pmdcl->GetClusSigmaX();
290 clusdata[5] = pmdcl->GetClusSigmaY();
291
292 AddRecPoint(idet,ismn,clusdata);
293
294 Int_t ncell = (Int_t) clusdata[3];
27a78deb 295 if (ncell > 19) ncell = 19;
8fbad6d3 296 for(Int_t ihit = 0; ihit < ncell; ihit++)
297 {
298 Int_t celldataX = pmdcl->GetClusCellX(ihit);
299 Int_t celldataY = pmdcl->GetClusCellY(ihit);
920e13db 300 Int_t celldataTr = pmdcl->GetClusCellTrack(ihit);
301 Int_t celldataPid = pmdcl->GetClusCellPid(ihit);
b0e4d1e1 302 Float_t celldataAdc = pmdcl->GetClusCellAdc(ihit);
303 AddRecHit(celldataX, celldataY, celldataTr, celldataPid, celldataAdc);
8fbad6d3 304 }
305 branch2->Fill();
306 ResetRechit();
307 }
26e08752 308 pmdcont->Delete();
920e13db 309
8fbad6d3 310 branch1->Fill();
311 ResetRecpoint();
312
313 } // modules
314
920e13db 315
8fbad6d3 316 ResetCellADC();
8fbad6d3 317
318 // delete the pointers
319 delete pmdclust;
320 delete pmdcont;
8fbad6d3 321}
322// ------------------------------------------------------------------------- //
323
7e9508a7 324void AliPMDClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
325 TTree *clustersTree)
326{
327 // Converts RAW data to recpoints after running clustering
328 // algorithm on CPV and PREshower plane
329 //
2c1131dd 330 // This method is called at the time of reconstruction from RAW data
2332574a 331
7e9508a7 332
920e13db 333 AliPMDddldata *pmdddl = 0x0;
334 AliPMDcluster *pmdcl = 0x0;
5254fcca 335
5c5cadd5 336 Float_t clusdata[6];
df42ab21 337 TObjArray pmdddlcont;
7e9508a7 338
339 TObjArray *pmdcont = new TObjArray();
5254fcca 340
5c5cadd5 341 AliPMDClustering *pmdclust = new AliPMDClusteringV1();
7e9508a7 342
08cf9354 343 // open the ddl file info to know the module
344 TString ddlinfofileName(gSystem->Getenv("ALICE_ROOT"));
345 ddlinfofileName += "/PMD/PMD_ddl_info.dat";
346
347 ifstream infileddl;
348 infileddl.open(ddlinfofileName.Data(), ios::in); // ascii file
349 if(!infileddl) AliError("Could not read the ddl info file");
350
be8b7039 351 Int_t ddlno = 0;
352 Int_t modno = 0;
353 Int_t modulePerDDL = 0;
354 Int_t moduleddl[6] = {0,0,0,0,0,0};
08cf9354 355
356 for(Int_t jddl = 0; jddl < 6; jddl++)
357 {
358 if (infileddl.eof()) break;
359 infileddl >> ddlno >> modulePerDDL;
360 moduleddl[jddl] = modulePerDDL;
361
362 if (modulePerDDL == 0) continue;
363 for (Int_t im = 0; im < modulePerDDL; im++)
364 {
365 infileddl >> modno;
366 }
367 }
368
369 infileddl.close();
370
b9e2f33c 371 // Set the minimum noise cut per module before clustering
372
373 fRecoParam = AliPMDReconstructor::GetRecoParam();
374
375 if(fRecoParam == 0x0)
376 {
377 AliFatal("No Reco Param found for PMD!!!");
378 }
379
7e9508a7 380 ResetRecpoint();
381
382 Int_t bufsize = 16000;
5c5cadd5 383 TBranch *branch1 = clustersTree->Branch("PMDRecpoint", &fRecpoints, bufsize);
384
385 TBranch * branch2 = clustersTree->Branch("PMDRechit", &fRechits, bufsize);
7e9508a7 386
7e9508a7 387 const Int_t kRow = 48;
388 const Int_t kCol = 96;
389
390 Int_t idet = 0;
391 Int_t iSMN = 0;
2332574a 392
722ccc67 393 Int_t indexDDL = -1;
394 AliPMDRawStream pmdinput(rawReader);
8fbad6d3 395
722ccc67 396 while ((indexDDL = pmdinput.DdlData(&pmdddlcont)) >=0)
08cf9354 397 {
398 iSMN = moduleddl[indexDDL];
399
7e9508a7 400 Int_t ***precpvADC;
401 precpvADC = new int **[iSMN];
402 for (Int_t i=0; i<iSMN; i++) precpvADC[i] = new int *[kRow];
403 for (Int_t i=0; i<iSMN;i++)
404 {
405 for (Int_t j=0; j<kRow; j++) precpvADC[i][j] = new int [kCol];
406 }
407 for (Int_t i = 0; i < iSMN; i++)
408 {
409 for (Int_t j = 0; j < kRow; j++)
410 {
411 for (Int_t k = 0; k < kCol; k++)
412 {
413 precpvADC[i][j][k] = 0;
414 }
415 }
416 }
417 ResetCellADC();
2332574a 418
df42ab21 419 Int_t indexsmn = 0;
420 Int_t ientries = pmdddlcont.GetEntries();
421 for (Int_t ient = 0; ient < ientries; ient++)
7e9508a7 422 {
920e13db 423 pmdddl = (AliPMDddldata*)pmdddlcont.UncheckedAt(ient);
7e9508a7 424
df42ab21 425 Int_t det = pmdddl->GetDetector();
426 Int_t smn = pmdddl->GetSMN();
427 //Int_t mcm = pmdddl->GetMCM();
428 //Int_t chno = pmdddl->GetChannel();
429 Int_t row = pmdddl->GetRow();
430 Int_t col = pmdddl->GetColumn();
431 Int_t sig = pmdddl->GetSignal();
2c1131dd 432
72f10132 433
8a523d22 434 if(det < 0 || det > 1)
435 {
436 AliError(Form("*CPV/PRE NUMBER WRONG %d *",det));
437 continue;
438 }
72f10132 439 if(smn < 0 || smn > 23)
2c1131dd 440 {
441 AliError(Form("*MODULE NUMBER WRONG %d *",smn));
442 continue;
443 }
72f10132 444 if(row < 0 || row > 47 || col < 0 || col > 95)
f0290c15 445 {
446 AliError(Form("*Row %d and Column NUMBER %d NOT Valid *",
447 row, col));
08cf9354 448
f0290c15 449 continue;
450 }
2c1131dd 451
72f10132 452
35535af7 453 // Pedestal Subtraction
454 Int_t pedmeanrms = fCalibPed->GetPedMeanRms(det,smn,row,col);
37edc588 455 Int_t pedrms1 = (Int_t) pedmeanrms%100;
35535af7 456 Float_t pedrms = (Float_t)pedrms1/10.;
457 Float_t pedmean = (Float_t) (pedmeanrms - pedrms1)/1000.0;
458
459 //printf("%f %f\n",pedmean, pedrms);
460
461 // Float_t sig1 = (Float_t) sig;
462 Float_t sig1 = (Float_t) sig - (pedmean + 3.0*pedrms);
463
8a523d22 464
7adbd8e9 465 // Hot cell - set the cell adc = 0
466 Float_t hotflag = fCalibHot->GetHotChannel(det,smn,row,col);
8a523d22 467 if (hotflag == 1.) sig1 = 0;
7adbd8e9 468
09a06455 469 // CALIBRATION
35535af7 470 Float_t gain = fCalibGain->GetGainFact(det,smn,row,col);
df42ab21 471 //printf("sig = %d gain = %f\n",sig,gain);
09a06455 472 sig = (Int_t) (sig1*gain);
7e9508a7 473
08cf9354 474 if (indexDDL == 0)
475 {
476 if (det != 0)
477 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
478 indexDDL, det));
479 if (iSMN == 6)
480 {
481 indexsmn = smn;
482 }
483 else if (iSMN == 12)
484 {
485 if (smn < 6)
486 indexsmn = smn;
487 else if (smn >= 18 && smn < 24)
488 indexsmn = smn-12;
489 }
490 }
491 else if (indexDDL >= 1 && indexDDL < 4)
7e9508a7 492 {
493 if (det != 0)
494 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
495 indexDDL, det));
496 indexsmn = smn - indexDDL * 6;
497 }
498 else if (indexDDL == 4)
499 {
500 if (det != 1)
501 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
502 indexDDL, det));
503 if (smn < 6)
504 {
505 indexsmn = smn;
506 }
df42ab21 507 else if (smn >= 18 && smn < 24)
7e9508a7 508 {
df42ab21 509 indexsmn = smn - 12;
7e9508a7 510 }
511 }
512 else if (indexDDL == 5)
513 {
514 if (det != 1)
515 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
516 indexDDL, det));
08cf9354 517 if (smn >= 6 && smn < 18)
7e9508a7 518 {
df42ab21 519 indexsmn = smn - 6;
7e9508a7 520 }
521 }
08cf9354 522
7e9508a7 523 precpvADC[indexsmn][row][col] = sig;
df42ab21 524 }
525
26e08752 526 pmdddlcont.Delete();
7e9508a7 527
5254fcca 528 Int_t totAdcMod = 0;
529
7e9508a7 530 Int_t ismn = 0;
78fc1b96 531 for (indexsmn = 0; indexsmn < iSMN; indexsmn++)
7e9508a7 532 {
533 ResetCellADC();
5254fcca 534 totAdcMod = 0;
7e9508a7 535 for (Int_t irow = 0; irow < kRow; irow++)
536 {
537 for (Int_t icol = 0; icol < kCol; icol++)
538 {
920e13db 539 fCellTrack[irow][icol] = -1;
540 fCellPid[irow][icol] = -1;
541
7e9508a7 542 fCellADC[irow][icol] =
543 (Double_t) precpvADC[indexsmn][irow][icol];
5254fcca 544 totAdcMod += precpvADC[indexsmn][irow][icol];
7e9508a7 545 } // row
546 } // col
df42ab21 547
08cf9354 548 if (indexDDL == 0)
549 {
550 if (iSMN == 6)
551 {
552 ismn = indexsmn;
553 }
554 else if (iSMN == 12)
555 {
556
557 if (indexsmn < 6)
558 ismn = indexsmn;
559 else if (indexsmn >= 6 && indexsmn < 12)
560 ismn = indexsmn + 12;
561 }
562 idet = 0;
563 }
564 else if (indexDDL >= 1 && indexDDL < 4)
7e9508a7 565 {
566 ismn = indexsmn + indexDDL * 6;
567 idet = 0;
568 }
569 else if (indexDDL == 4)
570 {
571 if (indexsmn < 6)
572 {
573 ismn = indexsmn;
574 }
575 else if (indexsmn >= 6 && indexsmn < 12)
576 {
df42ab21 577 ismn = indexsmn + 12;
7e9508a7 578 }
579 idet = 1;
580 }
581 else if (indexDDL == 5)
582 {
08cf9354 583 ismn = indexsmn + 6;
7e9508a7 584 idet = 1;
585 }
586
5254fcca 587 if (totAdcMod <= 0) continue;
920e13db 588
b9e2f33c 589 Int_t imod = idet*24 + ismn;
590
963fe679 591 // Int_t cluspar = fRecoParam->GetPbPbParam()->GetClusteringParam();
4755c3f1 592 Int_t cluspar = fRecoParam->GetPPParam()->GetClusteringParam();
963fe679 593 // Int_t cluspar = fRecoParam->GetCosmicParam()->GetClusteringParam();
4755c3f1 594 pmdclust->SetClusteringParam(cluspar);
963fe679 595 Float_t encut = fNoiseCut->GetNoiseCut(imod);
b9e2f33c 596
963fe679 597 pmdclust->SetEdepCut(encut);
27a78deb 598 pmdclust->DoClust(idet,ismn,fCellTrack,fCellPid,fCellADC,pmdcont);
920e13db 599
7e9508a7 600 Int_t nentries1 = pmdcont->GetEntries();
601
602 AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
603
604 for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
605 {
920e13db 606 pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
7e9508a7 607 idet = pmdcl->GetDetector();
608 ismn = pmdcl->GetSMN();
609 clusdata[0] = pmdcl->GetClusX();
610 clusdata[1] = pmdcl->GetClusY();
611 clusdata[2] = pmdcl->GetClusADC();
612 clusdata[3] = pmdcl->GetClusCells();
5c5cadd5 613 clusdata[4] = pmdcl->GetClusSigmaX();
614 clusdata[5] = pmdcl->GetClusSigmaY();
7e9508a7 615
616 AddRecPoint(idet,ismn,clusdata);
09a06455 617
5c5cadd5 618 Int_t ncell = (Int_t) clusdata[3];
27a78deb 619 if (ncell > 19) ncell = 19;
5c5cadd5 620 for(Int_t ihit = 0; ihit < ncell; ihit++)
621 {
622 Int_t celldataX = pmdcl->GetClusCellX(ihit);
623 Int_t celldataY = pmdcl->GetClusCellY(ihit);
920e13db 624 Int_t celldataTr = pmdcl->GetClusCellTrack(ihit);
625 Int_t celldataPid = pmdcl->GetClusCellPid(ihit);
b0e4d1e1 626 Float_t celldataAdc = pmdcl->GetClusCellAdc(ihit);
627 AddRecHit(celldataX, celldataY, celldataTr, celldataPid, celldataAdc);
5c5cadd5 628 }
629 branch2->Fill();
630 ResetRechit();
631
7e9508a7 632 }
26e08752 633 pmdcont->Delete();
920e13db 634
5c5cadd5 635 branch1->Fill();
7e9508a7 636 ResetRecpoint();
637
638
639 } // smn
640
641 for (Int_t i=0; i<iSMN; i++)
642 {
643 for (Int_t j=0; j<kRow; j++) delete [] precpvADC[i][j];
644 }
645 for (Int_t i=0; i<iSMN; i++) delete [] precpvADC[i];
c9e0bd24 646 delete [] precpvADC;
722ccc67 647
7e9508a7 648 } // DDL Loop
920e13db 649
7e9508a7 650
651 ResetCellADC();
652
653 // delete the pointers
654 delete pmdclust;
655 delete pmdcont;
01709453 656}
1758e4fe 657// ------------------------------------------------------------------------- //
1758e4fe 658void AliPMDClusterFinder::AddRecPoint(Int_t idet,Int_t ismn,Float_t *clusdata)
01709453 659{
b208c6a3 660 // Add Reconstructed points
661 //
01709453 662 TClonesArray &lrecpoints = *fRecpoints;
ed228cbc 663 AliPMDrecpoint1 *newrecpoint;
1758e4fe 664 newrecpoint = new AliPMDrecpoint1(idet, ismn, clusdata);
ed228cbc 665 new(lrecpoints[fNpoint++]) AliPMDrecpoint1(newrecpoint);
01709453 666 delete newrecpoint;
667}
1758e4fe 668// ------------------------------------------------------------------------- //
920e13db 669void AliPMDClusterFinder::AddRecHit(Int_t celldataX,Int_t celldataY,
b0e4d1e1 670 Int_t celldataTr, Int_t celldataPid,
671 Float_t celldataAdc)
5c5cadd5 672{
673 // Add associated cell hits to the Reconstructed points
674 //
675 TClonesArray &lrechits = *fRechits;
676 AliPMDrechit *newrechit;
b0e4d1e1 677 newrechit = new AliPMDrechit(celldataX, celldataY, celldataTr, celldataPid, celldataAdc);
5c5cadd5 678 new(lrechits[fNhit++]) AliPMDrechit(newrechit);
679 delete newrechit;
680}
681// ------------------------------------------------------------------------- //
01709453 682void AliPMDClusterFinder::ResetCellADC()
683{
b208c6a3 684 // Reset the individual cell ADC value to zero
685 //
5e6a9312 686 for(Int_t irow = 0; irow < fgkRow; irow++)
01709453 687 {
5e6a9312 688 for(Int_t icol = 0; icol < fgkCol; icol++)
01709453 689 {
920e13db 690 fCellTrack[irow][icol] = -1;
691 fCellPid[irow][icol] = -1;
692 fCellADC[irow][icol] = 0.;
01709453 693 }
694 }
695}
1758e4fe 696// ------------------------------------------------------------------------- //
01709453 697void AliPMDClusterFinder::ResetRecpoint()
698{
b208c6a3 699 // Clear the list of reconstructed points
01709453 700 fNpoint = 0;
701 if (fRecpoints) fRecpoints->Clear();
702}
1758e4fe 703// ------------------------------------------------------------------------- //
5c5cadd5 704void AliPMDClusterFinder::ResetRechit()
705{
706 // Clear the list of reconstructed points
707 fNhit = 0;
708 if (fRechits) fRechits->Clear();
709}
710// ------------------------------------------------------------------------- //
55601d47 711void AliPMDClusterFinder::Load()
712{
ebd83c56 713 // Load all the *.root files
55601d47 714 //
715 fPMDLoader->LoadDigits("READ");
716 fPMDLoader->LoadRecPoints("recreate");
717}
718// ------------------------------------------------------------------------- //
01b56f5c 719void AliPMDClusterFinder::LoadClusters()
720{
721 // Load all the *.root files
722 //
723 fPMDLoader->LoadRecPoints("recreate");
724}
725// ------------------------------------------------------------------------- //
ebd83c56 726void AliPMDClusterFinder::UnLoad()
01709453 727{
b208c6a3 728 // Unload all the *.root files
729 //
ebd83c56 730 fPMDLoader->UnloadDigits();
731 fPMDLoader->UnloadRecPoints();
01709453 732}
1758e4fe 733// ------------------------------------------------------------------------- //
01b56f5c 734void AliPMDClusterFinder::UnLoadClusters()
735{
736 // Unload all the *.root files
737 //
738 fPMDLoader->UnloadRecPoints();
739}
740// ------------------------------------------------------------------------- //
35535af7 741AliPMDCalibData* AliPMDClusterFinder::GetCalibGain() const
09a06455 742{
743 // The run number will be centralized in AliCDBManager,
744 // you don't need to set it here!
745 // Added by ZA
35535af7 746 AliCDBEntry *entry = AliCDBManager::Instance()->Get("PMD/Calib/Gain");
09a06455 747
0dd3d6f9 748 if(!entry) AliFatal("Calibration object retrieval failed! ");
09a06455 749
750 AliPMDCalibData *calibdata=0;
751 if (entry) calibdata = (AliPMDCalibData*) entry->GetObject();
752
0dd3d6f9 753 if (!calibdata) AliFatal("No calibration data from calibration database !");
09a06455 754
755 return calibdata;
756}
35535af7 757// ------------------------------------------------------------------------- //
35535af7 758AliPMDPedestal* AliPMDClusterFinder::GetCalibPed() const
759{
760 // The run number will be centralized in AliCDBManager,
761 // you don't need to set it here!
762 AliCDBEntry *entry = AliCDBManager::Instance()->Get("PMD/Calib/Ped");
763
764 if(!entry) AliFatal("Pedestal object retrieval failed!");
765
766 AliPMDPedestal *pedestal = 0;
767 if (entry) pedestal = (AliPMDPedestal*) entry->GetObject();
768
769 if (!pedestal) AliFatal("No pedestal data from pedestal database !");
770
771 return pedestal;
772}
7adbd8e9 773//--------------------------------------------------------------------//
774AliPMDHotData* AliPMDClusterFinder::GetCalibHot() const
775{
776 // The run number will be centralized in AliCDBManager,
777 // you don't need to set it here!
778 AliCDBEntry *entry = AliCDBManager::Instance()->Get("PMD/Calib/Hot");
779
780 if(!entry) AliFatal("HotData object retrieval failed!");
781
782 AliPMDHotData *hot = 0;
783 if (entry) hot = (AliPMDHotData*) entry->GetObject();
784
785 if (!hot) AliFatal("No hot data from database !");
786
787 return hot;
788}
963fe679 789//--------------------------------------------------------------------//
790AliPMDNoiseCut* AliPMDClusterFinder::GetNoiseCut() const
791{
792 // The run number will be centralized in AliCDBManager,
793 // you don't need to set it here!
794 AliCDBEntry *entry = AliCDBManager::Instance()->Get("PMD/Calib/NoiseCut");
795
796 if(!entry) AliFatal("Noisecut object retrieval failed!");
797
798 AliPMDNoiseCut *ncut = 0;
799 if (entry) ncut = (AliPMDNoiseCut*) entry->GetObject();
800
801 if (!ncut) AliFatal("No noise cut data from database !");
802
803 return ncut;
804}
7adbd8e9 805