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