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