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