]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDtracker.cxx
Bugfix
[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 "AliESD.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()
74 {
75   // Destructor
76   if (fRecpoints)
77     {
78       fRecpoints->Delete();
79       delete fRecpoints;
80       fRecpoints=0;
81     }
82   if (fPMDcontin)
83     {
84       fPMDcontin->Delete();
85       delete fPMDcontin;
86       fPMDcontin=0;
87     }
88   if (fPMDcontout)
89     {
90       fPMDcontout->Delete();
91       delete fPMDcontout;
92       fPMDcontout=0;
93     }
94 }
95 //--------------------------------------------------------------------//
96 void AliPMDtracker::LoadClusters(TTree *treein)
97 {
98   // Load the Reconstructed tree
99   fTreeR = treein;
100 }
101 //--------------------------------------------------------------------//
102 void AliPMDtracker::Clusters2Tracks(AliESD *event)
103 {
104   // Converts digits to recpoints after running clustering
105   // algorithm on CPV plane and PREshower plane
106   //
107
108   Int_t   idet;
109   Int_t   ismn;
110   Float_t clusdata[6];
111
112   TBranch *branch = fTreeR->GetBranch("PMDRecpoint");
113   if (!branch)
114     {
115       AliError("PMDRecpoint branch not found");
116       return;
117     }
118   branch->SetAddress(&fRecpoints);  
119   
120   Int_t   nmodules = (Int_t) branch->GetEntries();
121   
122   AliDebug(1,Form("Number of modules filled in treeR = %d",nmodules));
123   for (Int_t imodule = 0; imodule < nmodules; imodule++)
124     {
125       branch->GetEntry(imodule); 
126       Int_t nentries = fRecpoints->GetLast();
127       AliDebug(2,Form("Number of clusters per modules filled in treeR = %d"
128                       ,nentries));
129       for(Int_t ient = 0; ient < nentries+1; ient++)
130         {
131           fPMDrecpoint = (AliPMDrecpoint1*)fRecpoints->UncheckedAt(ient);
132           idet        = fPMDrecpoint->GetDetector();
133           ismn        = fPMDrecpoint->GetSMNumber();
134           clusdata[0] = fPMDrecpoint->GetClusX();
135           clusdata[1] = fPMDrecpoint->GetClusY();
136           clusdata[2] = fPMDrecpoint->GetClusADC();
137           clusdata[3] = fPMDrecpoint->GetClusCells();
138           clusdata[4] = fPMDrecpoint->GetClusSigmaX();
139           clusdata[5] = fPMDrecpoint->GetClusSigmaY();
140
141           fPMDclin = new AliPMDrecpoint1(idet,ismn,clusdata);
142           fPMDcontin->Add(fPMDclin);
143         }
144     }
145
146   AliPMDDiscriminator *pmddiscriminator = new AliPMDEmpDiscriminator();
147   pmddiscriminator->Discrimination(fPMDcontin,fPMDcontout);
148
149   const Float_t kzpos = 361.5;    // middle of the PMD
150   Int_t   ism =0, ium=0;
151   Int_t   det,smn;
152   Float_t xpos,ypos;
153   Float_t xpad = 0, ypad = 0;
154   Float_t adc, ncell, rad;
155   Float_t xglobal, yglobal, zglobal = 0;
156   Float_t pid;
157
158
159   Int_t nentries2 = fPMDcontout->GetEntries();
160   AliDebug(1,Form("Number of clusters coming after discrimination = %d"
161                   ,nentries2));
162   for (Int_t ient1 = 0; ient1 < nentries2; ient1++)
163     {
164       fPMDclout = (AliPMDclupid*)fPMDcontout->UncheckedAt(ient1);
165       
166       det   = fPMDclout->GetDetector();
167       smn   = fPMDclout->GetSMN();
168       xpos  = fPMDclout->GetClusX();
169       ypos  = fPMDclout->GetClusY();
170       adc   = fPMDclout->GetClusADC();
171       ncell = fPMDclout->GetClusCells();
172       rad   = fPMDclout->GetClusRadius();
173       pid   = fPMDclout->GetClusPID();
174       
175       //
176       // Now change the xpos and ypos to its original values
177       // for the unit modules which are earlier changed.
178       // xpad and ypad are the real positions.
179       //
180       /**********************************************************************
181        *    det   : Detector, 0: PRE & 1:CPV                                *
182        *    smn   : Serial Module Number from which Super Module Number     *
183        *            and Unit Module Numbers are extracted                   *
184        *    xpos  : x-position of the cluster                               *
185        *    ypos  : y-position of the cluster                               *
186        *            THESE xpos & ypos are not the true xpos and ypos        *
187        *            for some of the unit modules. They are rotated.         *
188        *    adc   : ADC contained in the cluster                            *
189        *    ncell : Number of cells contained in the cluster                *
190        *    rad   : radius of the cluster (1d fit)                          *
191        *    ism   : Supermodule number extracted from smn                   *
192        *    ium   : Unit module number extracted from smn                   *
193        *    xpad  : TRUE x-position of the cluster                          *
194        *    ypad  : TRUE y-position of the cluster                          *
195        **********************************************************************/
196       //
197       if(det == 0 || det == 1)
198         {
199           if(smn < 12)
200             {
201               ism  = smn/6;
202               ium  = smn - ism*6;
203               xpad = ypos;
204               ypad = xpos;
205             }
206           else if( smn >= 12 && smn < 24)
207             {
208               ism  = smn/6;
209               ium  = smn - ism*6;
210               xpad = xpos;
211               ypad = ypos;
212             }
213         }
214      
215       fPMDutil->RectGeomCellPos(ism,ium,xpad,ypad,xglobal,yglobal);
216
217       if (det == 0)
218         {
219           zglobal = kzpos + 1.6; // PREshower plane
220         }
221       else if (det == 1)
222         {
223           zglobal = kzpos - 1.7; // CPV plane
224         }
225
226
227       // Fill ESD
228
229       AliESDPmdTrack *esdpmdtr = new  AliESDPmdTrack();
230
231       esdpmdtr->SetDetector(det);
232       esdpmdtr->SetClusterX(xglobal);
233       esdpmdtr->SetClusterY(yglobal);
234       esdpmdtr->SetClusterZ(zglobal);
235       esdpmdtr->SetClusterADC(adc);
236       esdpmdtr->SetClusterCells(ncell);
237       esdpmdtr->SetClusterPID(pid);
238
239       event->AddPmdTrack(esdpmdtr);
240     }
241 }
242 //--------------------------------------------------------------------//
243 void AliPMDtracker::SetVertex(Double_t vtx[3], Double_t evtx[3])
244 {
245   fXvertex = vtx[0];
246   fYvertex = vtx[1];
247   fZvertex = vtx[2];
248   fSigmaX  = evtx[0];
249   fSigmaY  = evtx[1];
250   fSigmaZ  = evtx[2];
251 }
252 //--------------------------------------------------------------------//
253 void AliPMDtracker::ResetClusters()
254 {
255   if (fRecpoints) fRecpoints->Clear();
256 }
257 //--------------------------------------------------------------------//