]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDtracker.cxx
readers updated (mini header -> data header)
[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   if (!branch) return;
115   branch->SetAddress(&fRecpoints);  
116   
117   Int_t   nmodules = (Int_t) fTreeR->GetEntries();
118   cout << " nmodules = " << nmodules << endl;
119   for (Int_t imodule = 0; imodule < nmodules; imodule++)
120     {
121       fTreeR->GetEntry(imodule); 
122       Int_t nentries = fRecpoints->GetLast();
123       //      cout << " nentries = " << nentries << endl;
124       for(Int_t ient = 0; ient < nentries+1; ient++)
125         {
126           fPMDrecpoint = (AliPMDrecpoint1*)fRecpoints->UncheckedAt(ient);
127           idet        = fPMDrecpoint->GetDetector();
128           ismn        = fPMDrecpoint->GetSMNumber();
129           clusdata[0] = fPMDrecpoint->GetClusX();
130           clusdata[1] = fPMDrecpoint->GetClusY();
131           clusdata[2] = fPMDrecpoint->GetClusADC();
132           clusdata[3] = fPMDrecpoint->GetClusCells();
133           clusdata[4] = fPMDrecpoint->GetClusRadius();
134           
135           fPMDclin = new AliPMDcluster(idet,ismn,clusdata);
136           fPMDcontin->Add(fPMDclin);
137         }
138     }
139
140   fPMDdiscriminator->Discrimination(fPMDcontin,fPMDcontout);
141
142   const Float_t kzpos = 361.5;
143   Int_t   ism =0, ium=0;
144   Int_t   det,smn;
145   Float_t xpos,ypos;
146   Float_t xpad = 0, ypad = 0;
147   Float_t adc, ncell, rad;
148   Float_t xglobal, yglobal;
149   Float_t pid;
150
151   Float_t zglobal = kzpos + (Float_t) fZvertex;
152
153   Int_t nentries2 = fPMDcontout->GetEntries();
154   cout << " nentries2 = " << nentries2 << endl;
155   for (Int_t ient1 = 0; ient1 < nentries2; ient1++)
156     {
157       fPMDclout = (AliPMDclupid*)fPMDcontout->UncheckedAt(ient1);
158       
159       det   = fPMDclout->GetDetector();
160       smn   = fPMDclout->GetSMN();
161       xpos  = fPMDclout->GetClusX();
162       ypos  = fPMDclout->GetClusY();
163       adc   = fPMDclout->GetClusADC();
164       ncell = fPMDclout->GetClusCells();
165       rad   = fPMDclout->GetClusRadius();
166       pid   = fPMDclout->GetClusPID();
167       
168       //
169       // Now change the xpos and ypos to its original values
170       // for the unit modules which are earlier changed.
171       // xpad and ypad are the real positions.
172       //
173       /**********************************************************************
174        *    det   : Detector, 0: PRE & 1:CPV                                *
175        *    smn   : Serial Module Number from which Super Module Number     *
176        *            and Unit Module Numbers are extracted                   *
177        *    xpos  : x-position of the cluster                               *
178        *    ypos  : y-position of the cluster                               *
179        *            THESE xpos & ypos are not the true xpos and ypos        *
180        *            for some of the unit modules. They are rotated.         *
181        *    adc   : ADC contained in the cluster                            *
182        *    ncell : Number of cells contained in the cluster                *
183        *    rad   : radius of the cluster (1d fit)                          *
184        *    ism   : Supermodule number extracted from smn                   *
185        *    ium   : Unit module number extracted from smn                   *
186        *    xpad  : TRUE x-position of the cluster                          *
187        *    ypad  : TRUE y-position of the cluster                          *
188        **********************************************************************/
189       //
190       if(det == 0 || det == 1)
191         {
192           if(smn < 12)
193             {
194               ism  = smn/6;
195               ium  = smn - ism*6;
196               xpad = ypos;
197               ypad = xpos;
198             }
199           else if( smn >= 12 && smn < 24)
200             {
201               ism  = smn/6;
202               ium  = smn - ism*6;
203               xpad = xpos;
204               ypad = ypos;
205             }
206         }
207      
208       fPMDutil->RectGeomCellPos(ism,ium,xpad,ypad,xglobal,yglobal);
209       fPMDutil->SetXYZ(xglobal,yglobal,zglobal);
210       fPMDutil->CalculateEtaPhi();
211       Float_t theta = fPMDutil->GetTheta();
212       Float_t phi   = fPMDutil->GetPhi();
213
214       // Fill ESD
215
216       AliESDPmdTrack *esdpmdtr = new  AliESDPmdTrack();
217
218       esdpmdtr->SetDetector(det);
219       esdpmdtr->SetTheta(theta);
220       esdpmdtr->SetPhi(phi);
221       esdpmdtr->SetClusterADC(adc);
222       esdpmdtr->SetClusterPID(pid);
223
224       event->AddPmdTrack(esdpmdtr);
225     }
226
227 }
228 //--------------------------------------------------------------------//
229 void AliPMDtracker::SetVertex(Double_t vtx[3], Double_t evtx[3])
230 {
231   fXvertex = vtx[0];
232   fYvertex = vtx[1];
233   fZvertex = vtx[2];
234   fSigmaX  = evtx[0];
235   fSigmaY  = evtx[1];
236   fSigmaZ  = evtx[2];
237 }
238 //--------------------------------------------------------------------//
239 void AliPMDtracker::SetDebug(Int_t idebug)
240 {
241   fDebug = idebug;
242 }
243 //--------------------------------------------------------------------//
244 void AliPMDtracker::ResetClusters()
245 {
246   if (fRecpoints) fRecpoints->Clear();
247 }
248 //--------------------------------------------------------------------//