]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDtracker.cxx
isolated cell serach included in AliPMDClusteringV1 and kept in recpoints for offline...
[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 <TTree.h>
28 #include <TObjArray.h>
29 #include <TClonesArray.h>
30 #include <TFile.h>
31 #include <TBranch.h>
32 #include <TNtuple.h>
33 #include <TParticle.h>
34
35 #include "AliPMDcluster.h"
36 #include "AliPMDclupid.h"
37 #include "AliPMDrecpoint1.h"
38 #include "AliPMDUtility.h"
39 #include "AliPMDDiscriminator.h"
40 #include "AliPMDEmpDiscriminator.h"
41 #include "AliPMDtracker.h"
42
43 #include "AliESDPmdTrack.h"
44 #include "AliESDEvent.h"
45 #include "AliLog.h"
46
47 ClassImp(AliPMDtracker)
48
49 AliPMDtracker::AliPMDtracker():
50   fTreeR(0),
51   fRecpoints(new TClonesArray("AliPMDrecpoint1", 10)),
52   fPMDcontin(new TObjArray()),
53   fPMDcontout(new TObjArray()),
54   fPMDutil(new AliPMDUtility()),
55   fPMDrecpoint(0),
56   fPMDclin(0),
57   fPMDclout(0),
58   fXvertex(0.),
59   fYvertex(0.),
60   fZvertex(0.),
61   fSigmaX(0.),
62   fSigmaY(0.),
63   fSigmaZ(0.)
64 {
65   //
66   // Default Constructor
67   //
68 }
69 //--------------------------------------------------------------------//
70 AliPMDtracker:: AliPMDtracker(const AliPMDtracker & /* tracker */):
71   TObject(/* tracker */),
72   fTreeR(0),
73   fRecpoints(NULL),
74   fPMDcontin(NULL),
75   fPMDcontout(NULL),
76   fPMDutil(NULL),
77   fPMDrecpoint(0),
78   fPMDclin(0),
79   fPMDclout(0),
80   fXvertex(0.),
81   fYvertex(0.),
82   fZvertex(0.),
83   fSigmaX(0.),
84   fSigmaY(0.),
85   fSigmaZ(0.)
86 {
87   // copy constructor
88   AliError("Copy constructor not allowed");
89 }
90
91 //--------------------------------------------------------------------//
92 AliPMDtracker& AliPMDtracker::operator=(const AliPMDtracker & /* tracker */)
93 {
94  // assignment operator
95   AliError("Assignment operator not allowed");
96   return *this;
97 }
98
99 //--------------------------------------------------------------------//
100 AliPMDtracker::~AliPMDtracker()
101 {
102   // Destructor
103   if (fRecpoints)
104     {
105       fRecpoints->Clear();
106     }
107   if (fPMDcontin)
108     {
109       fPMDcontin->Delete();
110       delete fPMDcontin;
111       fPMDcontin=0;
112       
113     }
114   if (fPMDcontout)
115   {
116       fPMDcontout->Delete();
117       delete fPMDcontout;
118       fPMDcontout=0;
119
120     }
121   delete fPMDutil;
122 }
123 //--------------------------------------------------------------------//
124 void AliPMDtracker::LoadClusters(TTree *treein)
125 {
126   // Load the Reconstructed tree
127   fTreeR = treein;
128 }
129 //--------------------------------------------------------------------//
130 void AliPMDtracker::Clusters2Tracks(AliESDEvent *event)
131 {
132   // Converts digits to recpoints after running clustering
133   // algorithm on CPV plane and PREshower plane
134   //
135
136   Int_t   idet;
137   Int_t   ismn;
138   Float_t clusdata[6];
139
140   TBranch *branch = fTreeR->GetBranch("PMDRecpoint");
141   if (!branch)
142     {
143       AliError("PMDRecpoint branch not found");
144       return;
145     }
146   branch->SetAddress(&fRecpoints);  
147   
148   Int_t   nmodules = (Int_t) branch->GetEntries();
149   
150   AliDebug(1,Form("Number of modules filled in treeR = %d",nmodules));
151   for (Int_t imodule = 0; imodule < nmodules; imodule++)
152     {
153       branch->GetEntry(imodule); 
154       Int_t nentries = fRecpoints->GetLast();
155       AliDebug(2,Form("Number of clusters per modules filled in treeR = %d"
156                       ,nentries));
157       for(Int_t ient = 0; ient < nentries+1; ient++)
158         {
159           fPMDrecpoint = (AliPMDrecpoint1*)fRecpoints->UncheckedAt(ient);
160           idet        = fPMDrecpoint->GetDetector();
161           ismn        = fPMDrecpoint->GetSMNumber();
162           clusdata[0] = fPMDrecpoint->GetClusX();
163           clusdata[1] = fPMDrecpoint->GetClusY();
164           clusdata[2] = fPMDrecpoint->GetClusADC();
165           clusdata[3] = fPMDrecpoint->GetClusCells();
166           clusdata[4] = fPMDrecpoint->GetClusSigmaX();
167           clusdata[5] = fPMDrecpoint->GetClusSigmaY();
168           
169           if (clusdata[4] != -99. && clusdata[5] != -99.)
170             { 
171               fPMDclin = new AliPMDrecpoint1(idet,ismn,clusdata);
172               fPMDcontin->Add(fPMDclin);
173             }
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   fPMDcontin->Delete();
246   fPMDcontout->Delete();
247
248 }
249 //--------------------------------------------------------------------//
250 void AliPMDtracker::SetVertex(Double_t vtx[3], Double_t evtx[3])
251 {
252   fXvertex = vtx[0];
253   fYvertex = vtx[1];
254   fZvertex = vtx[2];
255   fSigmaX  = evtx[0];
256   fSigmaY  = evtx[1];
257   fSigmaZ  = evtx[2];
258 }
259 //--------------------------------------------------------------------//
260 void AliPMDtracker::ResetClusters()
261 {
262   if (fRecpoints) fRecpoints->Clear();
263 }
264 //--------------------------------------------------------------------//