1 /***************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 //-----------------------------------------------------//
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) //
23 //-----------------------------------------------------//
25 #include <Riostream.h>
27 #include <TObjArray.h>
28 #include <TClonesArray.h>
31 #include "AliRunLoader.h"
32 #include "AliLoader.h"
33 #include "AliRawReader.h"
35 #include "AliPMDdigit.h"
36 #include "AliPMDClusterFinder.h"
37 #include "AliPMDClustering.h"
38 #include "AliPMDClusteringV1.h"
39 #include "AliPMDcluster.h"
40 #include "AliPMDrecpoint1.h"
41 #include "AliPMDrechit.h"
42 #include "AliPMDRawStream.h"
46 ClassImp(AliPMDClusterFinder)
48 AliPMDClusterFinder::AliPMDClusterFinder():
53 fDigits(new TClonesArray("AliPMDdigit", 1000)),
54 fRecpoints(new TClonesArray("AliPMDrecpoint1", 1000)),
55 fRechits(new TClonesArray("AliPMDrechit", 1000)),
64 // ------------------------------------------------------------------------- //
65 AliPMDClusterFinder::AliPMDClusterFinder(AliRunLoader* runLoader):
66 fRunLoader(runLoader),
67 fPMDLoader(runLoader->GetLoader("PMDLoader")),
70 fDigits(new TClonesArray("AliPMDdigit", 1000)),
71 fRecpoints(new TClonesArray("AliPMDrecpoint1", 1000)),
72 fRechits(new TClonesArray("AliPMDrechit", 1000)),
81 // ------------------------------------------------------------------------- //
82 AliPMDClusterFinder::~AliPMDClusterFinder()
104 // ------------------------------------------------------------------------- //
106 void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt)
108 // Converts digits to recpoints after running clustering
109 // algorithm on CPV plane and PREshower plane
111 Int_t det = 0,smn = 0;
118 TObjArray *pmdcont = new TObjArray();
119 AliPMDClustering *pmdclust = new AliPMDClusteringV1();
121 pmdclust->SetEdepCut(fEcut);
123 fRunLoader->GetEvent(ievt);
126 fTreeD = fPMDLoader->TreeD();
129 AliFatal("AliPMDClusterFinder: Can not get TreeD");
132 AliPMDdigit *pmddigit;
133 TBranch *branch = fTreeD->GetBranch("PMDDigit");
134 branch->SetAddress(&fDigits);
138 fTreeR = fPMDLoader->TreeR();
141 fPMDLoader->MakeTree("R");
142 fTreeR = fPMDLoader->TreeR();
145 Int_t bufsize = 16000;
146 TBranch * branch1 = fTreeR->Branch("PMDRecpoint", &fRecpoints, bufsize);
147 TBranch * branch2 = fTreeR->Branch("PMDRechit", &fRechits, bufsize);
149 Int_t nmodules = (Int_t) fTreeD->GetEntries();
151 for (Int_t imodule = 0; imodule < nmodules; imodule++)
154 fTreeD->GetEntry(imodule);
155 Int_t nentries = fDigits->GetLast();
156 for (Int_t ient = 0; ient < nentries+1; ient++)
158 pmddigit = (AliPMDdigit*)fDigits->UncheckedAt(ient);
160 det = pmddigit->GetDetector();
161 smn = pmddigit->GetSMNumber();
162 xpos = pmddigit->GetRow();
163 ypos = pmddigit->GetColumn();
164 adc = pmddigit->GetADC();
165 //Int_t trno = pmddigit->GetTrackNumber();
166 fCellADC[xpos][ypos] = (Double_t) adc;
171 pmdclust->DoClust(idet,ismn,fCellADC,pmdcont);
173 Int_t nentries1 = pmdcont->GetEntries();
175 AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
177 for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
179 AliPMDcluster *pmdcl = (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
180 idet = pmdcl->GetDetector();
181 ismn = pmdcl->GetSMN();
182 clusdata[0] = pmdcl->GetClusX();
183 clusdata[1] = pmdcl->GetClusY();
184 clusdata[2] = pmdcl->GetClusADC();
185 clusdata[3] = pmdcl->GetClusCells();
186 clusdata[4] = pmdcl->GetClusSigmaX();
187 clusdata[5] = pmdcl->GetClusSigmaY();
189 AddRecPoint(idet,ismn,clusdata);
191 // This part will be worked out
193 Int_t ncell = (Int_t) clusdata[3];
194 for(Int_t ihit = 0; ihit < ncell; ihit++)
196 Int_t celldataX = pmdcl->GetClusCellX(ihit);
197 Int_t celldataY = pmdcl->GetClusCellY(ihit);
198 AddRecHit(celldataX, celldataY);
211 fPMDLoader = fRunLoader->GetLoader("PMDLoader");
212 fPMDLoader->WriteRecPoints("OVERWRITE");
214 // delete the pointers
219 // ------------------------------------------------------------------------- //
221 void AliPMDClusterFinder::Digits2RecPoints(AliRawReader *rawReader,
224 // Converts RAW data to recpoints after running clustering
225 // algorithm on CPV and PREshower plane
230 TObjArray *pmdcont = new TObjArray();
231 AliPMDClustering *pmdclust = new AliPMDClusteringV1();
233 pmdclust->SetEdepCut(fEcut);
237 Int_t bufsize = 16000;
238 TBranch *branch1 = clustersTree->Branch("PMDRecpoint", &fRecpoints, bufsize);
240 TBranch * branch2 = clustersTree->Branch("PMDRechit", &fRechits, bufsize);
242 const Int_t kDDL = 6;
243 const Int_t kRow = 48;
244 const Int_t kCol = 96;
249 for (Int_t indexDDL = 0; indexDDL < kDDL; indexDDL++)
255 else if (indexDDL >= 4)
260 precpvADC = new int **[iSMN];
261 for (Int_t i=0; i<iSMN; i++) precpvADC[i] = new int *[kRow];
262 for (Int_t i=0; i<iSMN;i++)
264 for (Int_t j=0; j<kRow; j++) precpvADC[i][j] = new int [kCol];
266 for (Int_t i = 0; i < iSMN; i++)
268 for (Int_t j = 0; j < kRow; j++)
270 for (Int_t k = 0; k < kCol; k++)
272 precpvADC[i][j][k] = 0;
278 AliPMDRawStream pmdinput(rawReader);
279 rawReader->Select(12, indexDDL, indexDDL);
280 while(pmdinput.Next())
282 Int_t det = pmdinput.GetDetector();
283 Int_t smn = pmdinput.GetSMN();
284 //Int_t mcm = pmdinput.GetMCM();
285 //Int_t chno = pmdinput.GetChannel();
286 Int_t row = pmdinput.GetRow();
287 Int_t col = pmdinput.GetColumn();
288 Int_t sig = pmdinput.GetSignal();
295 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
297 indexsmn = smn - indexDDL * 6;
299 else if (indexDDL == 4)
302 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
308 else if (smn >= 12 && smn < 18)
313 else if (indexDDL == 5)
316 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
318 if (smn >= 6 && smn < 12)
322 else if (smn >= 18 && smn < 24)
327 precpvADC[indexsmn][row][col] = sig;
331 for (Int_t indexsmn = 0; indexsmn < iSMN; indexsmn++)
334 for (Int_t irow = 0; irow < kRow; irow++)
336 for (Int_t icol = 0; icol < kCol; icol++)
338 fCellADC[irow][icol] =
339 (Double_t) precpvADC[indexsmn][irow][icol];
344 ismn = indexsmn + indexDDL * 6;
347 else if (indexDDL == 4)
353 else if (indexsmn >= 6 && indexsmn < 12)
359 else if (indexDDL == 5)
365 else if (indexsmn >= 6 && indexsmn < 12)
367 ismn = indexsmn + 12;
373 pmdclust->DoClust(idet,ismn,fCellADC,pmdcont);
374 Int_t nentries1 = pmdcont->GetEntries();
376 AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
378 for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
380 AliPMDcluster *pmdcl =
381 (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
382 idet = pmdcl->GetDetector();
383 ismn = pmdcl->GetSMN();
384 clusdata[0] = pmdcl->GetClusX();
385 clusdata[1] = pmdcl->GetClusY();
386 clusdata[2] = pmdcl->GetClusADC();
387 clusdata[3] = pmdcl->GetClusCells();
388 clusdata[4] = pmdcl->GetClusSigmaX();
389 clusdata[5] = pmdcl->GetClusSigmaY();
391 AddRecPoint(idet,ismn,clusdata);
393 Int_t ncell = (Int_t) clusdata[3];
394 for(Int_t ihit = 0; ihit < ncell; ihit++)
396 Int_t celldataX = pmdcl->GetClusCellX(ihit);
397 Int_t celldataY = pmdcl->GetClusCellY(ihit);
398 AddRecHit(celldataX, celldataY);
412 for (Int_t i=0; i<iSMN; i++)
414 for (Int_t j=0; j<kRow; j++) delete [] precpvADC[i][j];
416 for (Int_t i=0; i<iSMN; i++) delete [] precpvADC[i];
422 // delete the pointers
427 // ------------------------------------------------------------------------- //
429 void AliPMDClusterFinder::Digits2RecPoints(Int_t ievt, AliRawReader *rawReader)
431 // Converts RAW data to recpoints after running clustering
432 // algorithm on CPV and PREshower plane
437 TObjArray *pmdcont = new TObjArray();
439 AliPMDClustering *pmdclust = new AliPMDClusteringV1();
441 pmdclust->SetEdepCut(fEcut);
443 fRunLoader->GetEvent(ievt);
447 fTreeR = fPMDLoader->TreeR();
450 fPMDLoader->MakeTree("R");
451 fTreeR = fPMDLoader->TreeR();
453 Int_t bufsize = 16000;
454 TBranch *branch1 = fTreeR->Branch("PMDRecpoint", &fRecpoints, bufsize);
455 TBranch *branch2 = fTreeR->Branch("PMDRechit", &fRechits, bufsize);
457 const Int_t kDDL = 6;
458 const Int_t kRow = 48;
459 const Int_t kCol = 96;
464 for (Int_t indexDDL = 0; indexDDL < kDDL; indexDDL++)
470 else if (indexDDL >= 4)
475 precpvADC = new int **[iSMN];
476 for (Int_t i=0; i<iSMN; i++) precpvADC[i] = new int *[kRow];
477 for (Int_t i=0; i<iSMN;i++)
479 for (Int_t j=0; j<kRow; j++) precpvADC[i][j] = new int [kCol];
481 for (Int_t i = 0; i < iSMN; i++)
483 for (Int_t j = 0; j < kRow; j++)
485 for (Int_t k = 0; k < kCol; k++)
487 precpvADC[i][j][k] = 0;
493 AliPMDRawStream pmdinput(rawReader);
494 rawReader->Select(12, indexDDL, indexDDL);
495 while(pmdinput.Next())
497 Int_t det = pmdinput.GetDetector();
498 Int_t smn = pmdinput.GetSMN();
499 //Int_t mcm = pmdinput.GetMCM();
500 //Int_t chno = pmdinput.GetChannel();
501 Int_t row = pmdinput.GetRow();
502 Int_t col = pmdinput.GetColumn();
503 Int_t sig = pmdinput.GetSignal();
510 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
512 indexsmn = smn - indexDDL * 6;
514 else if (indexDDL == 4)
517 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
523 else if (smn >= 12 && smn < 18)
528 else if (indexDDL == 5)
531 AliError(Form("*DDL %d and Detector NUMBER %d NOT MATCHING *",
533 if (smn >= 6 && smn < 12)
537 else if (smn >= 18 && smn < 24)
542 precpvADC[indexsmn][row][col] = sig;
546 for (Int_t indexsmn = 0; indexsmn < iSMN; indexsmn++)
549 for (Int_t irow = 0; irow < kRow; irow++)
551 for (Int_t icol = 0; icol < kCol; icol++)
553 fCellADC[irow][icol] =
554 (Double_t) precpvADC[indexsmn][irow][icol];
559 ismn = indexsmn + indexDDL * 6;
562 else if (indexDDL == 4)
568 else if (indexsmn >= 6 && indexsmn < 12)
574 else if (indexDDL == 5)
580 else if (indexsmn >= 6 && indexsmn < 12)
582 ismn = indexsmn + 12;
588 pmdclust->DoClust(idet,ismn,fCellADC,pmdcont);
589 Int_t nentries1 = pmdcont->GetEntries();
591 AliDebug(1,Form("Total number of clusters/module = %d",nentries1));
593 for (Int_t ient1 = 0; ient1 < nentries1; ient1++)
595 AliPMDcluster *pmdcl =
596 (AliPMDcluster*)pmdcont->UncheckedAt(ient1);
597 idet = pmdcl->GetDetector();
598 ismn = pmdcl->GetSMN();
599 clusdata[0] = pmdcl->GetClusX();
600 clusdata[1] = pmdcl->GetClusY();
601 clusdata[2] = pmdcl->GetClusADC();
602 clusdata[3] = pmdcl->GetClusCells();
603 clusdata[4] = pmdcl->GetClusSigmaX();
604 clusdata[5] = pmdcl->GetClusSigmaY();
606 AddRecPoint(idet,ismn,clusdata);
609 Int_t ncell = (Int_t) clusdata[3];
610 for(Int_t ihit = 0; ihit < ncell; ihit++)
612 Int_t celldataX = pmdcl->GetClusCellX(ihit);
613 Int_t celldataY = pmdcl->GetClusCellY(ihit);
614 AddRecHit(celldataX, celldataY);
628 for (Int_t i=0; i<iSMN; i++)
630 for (Int_t j=0; j<kRow; j++) delete [] precpvADC[i][j];
632 for (Int_t i=0; i<iSMN; i++) delete [] precpvADC[i];
638 fPMDLoader = fRunLoader->GetLoader("PMDLoader");
639 fPMDLoader->WriteRecPoints("OVERWRITE");
641 // delete the pointers
646 // ------------------------------------------------------------------------- //
647 void AliPMDClusterFinder::SetCellEdepCut(Float_t ecut)
651 // ------------------------------------------------------------------------- //
652 void AliPMDClusterFinder::AddRecPoint(Int_t idet,Int_t ismn,Float_t *clusdata)
654 // Add Reconstructed points
656 TClonesArray &lrecpoints = *fRecpoints;
657 AliPMDrecpoint1 *newrecpoint;
658 newrecpoint = new AliPMDrecpoint1(idet, ismn, clusdata);
659 new(lrecpoints[fNpoint++]) AliPMDrecpoint1(newrecpoint);
662 // ------------------------------------------------------------------------- //
663 void AliPMDClusterFinder::AddRecHit(Int_t celldataX,Int_t celldataY)
665 // Add associated cell hits to the Reconstructed points
667 TClonesArray &lrechits = *fRechits;
668 AliPMDrechit *newrechit;
669 newrechit = new AliPMDrechit(celldataX, celldataY);
670 new(lrechits[fNhit++]) AliPMDrechit(newrechit);
673 // ------------------------------------------------------------------------- //
674 void AliPMDClusterFinder::ResetCellADC()
676 // Reset the individual cell ADC value to zero
678 for(Int_t irow = 0; irow < fgkRow; irow++)
680 for(Int_t icol = 0; icol < fgkCol; icol++)
682 fCellADC[irow][icol] = 0.;
686 // ------------------------------------------------------------------------- //
688 void AliPMDClusterFinder::ResetRecpoint()
690 // Clear the list of reconstructed points
692 if (fRecpoints) fRecpoints->Clear();
694 // ------------------------------------------------------------------------- //
695 void AliPMDClusterFinder::ResetRechit()
697 // Clear the list of reconstructed points
699 if (fRechits) fRechits->Clear();
701 // ------------------------------------------------------------------------- //
702 void AliPMDClusterFinder::Load()
704 // Load all the *.root files
706 fPMDLoader->LoadDigits("READ");
707 fPMDLoader->LoadRecPoints("recreate");
709 // ------------------------------------------------------------------------- //
710 void AliPMDClusterFinder::LoadClusters()
712 // Load all the *.root files
714 fPMDLoader->LoadRecPoints("recreate");
716 // ------------------------------------------------------------------------- //
717 void AliPMDClusterFinder::UnLoad()
719 // Unload all the *.root files
721 fPMDLoader->UnloadDigits();
722 fPMDLoader->UnloadRecPoints();
724 // ------------------------------------------------------------------------- //
725 void AliPMDClusterFinder::UnLoadClusters()
727 // Unload all the *.root files
729 fPMDLoader->UnloadRecPoints();
731 // ------------------------------------------------------------------------- //