]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDtracker.cxx
Fixes for some mem-leaks: most changes where pretty basic (i.e. adding deletes).
[u/mrichter/AliRoot.git] / PMD / AliPMDtracker.cxx
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
16 //-----------------------------------------------------//
17 //                                                     //
18 //           Date   : March 25 2004                    //
19 //  This reads the file PMD.RecPoints.root(TreeR),     //
20 //  calls the Clustering algorithm and stores the      //
21 //  clustering output in PMD.RecPoints.root(TreeR)     // 
22 //                                                     //
23 //-----------------------------------------------------//
24
25 #include <Riostream.h>
26 #include <TMath.h>
27 #include <TBRIK.h>
28 #include <TNode.h>
29 #include <TTree.h>
30 #include <TGeometry.h>
31 #include <TObjArray.h>
32 #include <TClonesArray.h>
33 #include <TFile.h>
34 #include <TBranch.h>
35 #include <TNtuple.h>
36 #include <TParticle.h>
37
38 #include "AliPMDcluster.h"
39 #include "AliPMDclupid.h"
40 #include "AliPMDrecpoint1.h"
41 #include "AliPMDUtility.h"
42 #include "AliPMDDiscriminator.h"
43 #include "AliPMDEmpDiscriminator.h"
44 #include "AliPMDtracker.h"
45
46 #include "AliESDPmdTrack.h"
47 #include "AliESDEvent.h"
48 #include "AliLog.h"
49
50 ClassImp(AliPMDtracker)
51
52 AliPMDtracker::AliPMDtracker():
53   fTreeR(0),
54   fRecpoints(new TClonesArray("AliPMDrecpoint1", 1000)),
55   fPMDcontin(new TObjArray()),
56   fPMDcontout(new TObjArray()),
57   fPMDutil(new AliPMDUtility()),
58   fPMDrecpoint(0),
59   fPMDclin(0),
60   fPMDclout(0),
61   fXvertex(0.),
62   fYvertex(0.),
63   fZvertex(0.),
64   fSigmaX(0.),
65   fSigmaY(0.),
66   fSigmaZ(0.)
67 {
68   //
69   // Default Constructor
70   //
71 }
72 //--------------------------------------------------------------------//
73 AliPMDtracker:: AliPMDtracker(const AliPMDtracker & /* tracker */):
74   TObject(/* tracker */),
75   fTreeR(0),
76   fRecpoints(NULL),
77   fPMDcontin(NULL),
78   fPMDcontout(NULL),
79   fPMDutil(NULL),
80   fPMDrecpoint(0),
81   fPMDclin(0),
82   fPMDclout(0),
83   fXvertex(0.),
84   fYvertex(0.),
85   fZvertex(0.),
86   fSigmaX(0.),
87   fSigmaY(0.),
88   fSigmaZ(0.)
89 {
90   // copy constructor
91   AliError("Copy constructor not allowed");
92 }
93
94 //--------------------------------------------------------------------//
95 AliPMDtracker& AliPMDtracker::operator=(const AliPMDtracker & /* tracker */)
96 {
97  // assignment operator
98   AliError("Assignment operator not allowed");
99   return *this;
100 }
101
102 //--------------------------------------------------------------------//
103 AliPMDtracker::~AliPMDtracker()
104 {
105   // Destructor
106   if (fRecpoints)
107     {
108       fRecpoints->Delete();
109       delete fRecpoints;
110       fRecpoints=0;
111     }
112   if (fPMDcontin)
113     {
114       fPMDcontin->Delete();
115       delete fPMDcontin;
116       fPMDcontin=0;
117     }
118   if (fPMDcontout)
119     {
120       fPMDcontout->Delete();
121       delete fPMDcontout;
122       fPMDcontout=0;
123     }
124   delete fPMDutil;
125 }
126 //--------------------------------------------------------------------//
127 void AliPMDtracker::LoadClusters(TTree *treein)
128 {
129   // Load the Reconstructed tree
130   fTreeR = treein;
131 }
132 //--------------------------------------------------------------------//
133 void AliPMDtracker::Clusters2Tracks(AliESDEvent *event)
134 {
135   // Converts digits to recpoints after running clustering
136   // algorithm on CPV plane and PREshower plane
137   //
138
139   Int_t   idet;
140   Int_t   ismn;
141   Float_t clusdata[6];
142
143   TBranch *branch = fTreeR->GetBranch("PMDRecpoint");
144   if (!branch)
145     {
146       AliError("PMDRecpoint branch not found");
147       return;
148     }
149   branch->SetAddress(&fRecpoints);  
150   
151   Int_t   nmodules = (Int_t) branch->GetEntries();
152   
153   AliDebug(1,Form("Number of modules filled in treeR = %d",nmodules));
154   for (Int_t imodule = 0; imodule < nmodules; imodule++)
155     {
156       branch->GetEntry(imodule); 
157       Int_t nentries = fRecpoints->GetLast();
158       AliDebug(2,Form("Number of clusters per modules filled in treeR = %d"
159                       ,nentries));
160       for(Int_t ient = 0; ient < nentries+1; ient++)
161         {
162           fPMDrecpoint = (AliPMDrecpoint1*)fRecpoints->UncheckedAt(ient);
163           idet        = fPMDrecpoint->GetDetector();
164           ismn        = fPMDrecpoint->GetSMNumber();
165           clusdata[0] = fPMDrecpoint->GetClusX();
166           clusdata[1] = fPMDrecpoint->GetClusY();
167           clusdata[2] = fPMDrecpoint->GetClusADC();
168           clusdata[3] = fPMDrecpoint->GetClusCells();
169           clusdata[4] = fPMDrecpoint->GetClusSigmaX();
170           clusdata[5] = fPMDrecpoint->GetClusSigmaY();
171
172           fPMDclin = new AliPMDrecpoint1(idet,ismn,clusdata);
173           fPMDcontin->Add(fPMDclin);
174         }
175     }
176
177   AliPMDDiscriminator *pmddiscriminator = new AliPMDEmpDiscriminator();
178   pmddiscriminator->Discrimination(fPMDcontin,fPMDcontout);
179
180   const Float_t kzpos = 361.5;    // middle of the PMD
181
182   Int_t   det,smn;
183   Float_t xpos,ypos;
184   Float_t adc, ncell, rad;
185   Float_t xglobal = 0., yglobal = 0., zglobal = 0;
186   Float_t pid;
187
188
189   Int_t nentries2 = fPMDcontout->GetEntries();
190   AliDebug(1,Form("Number of clusters coming after discrimination = %d"
191                   ,nentries2));
192   for (Int_t ient1 = 0; ient1 < nentries2; ient1++)
193     {
194       fPMDclout = (AliPMDclupid*)fPMDcontout->UncheckedAt(ient1);
195       
196       det   = fPMDclout->GetDetector();
197       smn   = fPMDclout->GetSMN();
198       xpos  = fPMDclout->GetClusX();
199       ypos  = fPMDclout->GetClusY();
200       adc   = fPMDclout->GetClusADC();
201       ncell = fPMDclout->GetClusCells();
202       rad   = fPMDclout->GetClusRadius();
203       pid   = fPMDclout->GetClusPID();
204       
205       //
206       /**********************************************************************
207        *    det   : Detector, 0: PRE & 1:CPV                                *
208        *    smn   : Serial Module Number 0 to 23 for each plane             *
209        *    xpos  : x-position of the cluster                               *
210        *    ypos  : y-position of the cluster                               *
211        *            THESE xpos & ypos are not the true xpos and ypos        *
212        *            for some of the unit modules. They are rotated.         *
213        *    adc   : ADC contained in the cluster                            *
214        *    ncell : Number of cells contained in the cluster                *
215        *    rad   : radius of the cluster (1d fit)                          *
216        **********************************************************************/
217       //
218
219       fPMDutil->RectGeomCellPos(smn,xpos,ypos,xglobal,yglobal);
220
221       if (det == 0)
222         {
223           zglobal = kzpos + 1.6; // PREshower plane
224         }
225       else if (det == 1)
226         {
227           zglobal = kzpos - 1.7; // CPV plane
228         }
229
230       // Fill ESD
231
232       AliESDPmdTrack *esdpmdtr = new  AliESDPmdTrack();
233
234       esdpmdtr->SetDetector(det);
235       esdpmdtr->SetClusterX(xglobal);
236       esdpmdtr->SetClusterY(yglobal);
237       esdpmdtr->SetClusterZ(zglobal);
238       esdpmdtr->SetClusterADC(adc);
239       esdpmdtr->SetClusterCells(ncell);
240       esdpmdtr->SetClusterPID(pid);
241
242       event->AddPmdTrack(esdpmdtr);
243     }
244 }
245 //--------------------------------------------------------------------//
246 void AliPMDtracker::SetVertex(Double_t vtx[3], Double_t evtx[3])
247 {
248   fXvertex = vtx[0];
249   fYvertex = vtx[1];
250   fZvertex = vtx[2];
251   fSigmaX  = evtx[0];
252   fSigmaY  = evtx[1];
253   fSigmaZ  = evtx[2];
254 }
255 //--------------------------------------------------------------------//
256 void AliPMDtracker::ResetClusters()
257 {
258   if (fRecpoints) fRecpoints->Clear();
259 }
260 //--------------------------------------------------------------------//