]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PMD/AliPMDtracker.cxx
dimensions are changed because of different modules
[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
151   Int_t   det,smn;
152   Float_t xpos,ypos;
153   Float_t adc, ncell, rad;
154   Float_t xglobal = 0., yglobal = 0., zglobal = 0;
155   Float_t pid;
156
157
158   Int_t nentries2 = fPMDcontout->GetEntries();
159   AliDebug(1,Form("Number of clusters coming after discrimination = %d"
160                   ,nentries2));
161   for (Int_t ient1 = 0; ient1 < nentries2; ient1++)
162     {
163       fPMDclout = (AliPMDclupid*)fPMDcontout->UncheckedAt(ient1);
164       
165       det   = fPMDclout->GetDetector();
166       smn   = fPMDclout->GetSMN();
167       xpos  = fPMDclout->GetClusX();
168       ypos  = fPMDclout->GetClusY();
169       adc   = fPMDclout->GetClusADC();
170       ncell = fPMDclout->GetClusCells();
171       rad   = fPMDclout->GetClusRadius();
172       pid   = fPMDclout->GetClusPID();
173       
174       //
175       /**********************************************************************
176        *    det   : Detector, 0: PRE & 1:CPV                                *
177        *    smn   : Serial Module Number 0 to 23 for each plane             *
178        *    xpos  : x-position of the cluster                               *
179        *    ypos  : y-position of the cluster                               *
180        *            THESE xpos & ypos are not the true xpos and ypos        *
181        *            for some of the unit modules. They are rotated.         *
182        *    adc   : ADC contained in the cluster                            *
183        *    ncell : Number of cells contained in the cluster                *
184        *    rad   : radius of the cluster (1d fit)                          *
185        **********************************************************************/
186       //
187
188       fPMDutil->RectGeomCellPos(smn,xpos,ypos,xglobal,yglobal);
189
190       if (det == 0)
191         {
192           zglobal = kzpos + 1.6; // PREshower plane
193         }
194       else if (det == 1)
195         {
196           zglobal = kzpos - 1.7; // CPV plane
197         }
198
199       // Fill ESD
200
201       AliESDPmdTrack *esdpmdtr = new  AliESDPmdTrack();
202
203       esdpmdtr->SetDetector(det);
204       esdpmdtr->SetClusterX(xglobal);
205       esdpmdtr->SetClusterY(yglobal);
206       esdpmdtr->SetClusterZ(zglobal);
207       esdpmdtr->SetClusterADC(adc);
208       esdpmdtr->SetClusterCells(ncell);
209       esdpmdtr->SetClusterPID(pid);
210
211       event->AddPmdTrack(esdpmdtr);
212     }
213 }
214 //--------------------------------------------------------------------//
215 void AliPMDtracker::SetVertex(Double_t vtx[3], Double_t evtx[3])
216 {
217   fXvertex = vtx[0];
218   fYvertex = vtx[1];
219   fZvertex = vtx[2];
220   fSigmaX  = evtx[0];
221   fSigmaY  = evtx[1];
222   fSigmaZ  = evtx[2];
223 }
224 //--------------------------------------------------------------------//
225 void AliPMDtracker::ResetClusters()
226 {
227   if (fRecpoints) fRecpoints->Clear();
228 }
229 //--------------------------------------------------------------------//