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