1 /***************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 //-----------------------------------------------------//
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) //
23 //-----------------------------------------------------//
25 #include <Riostream.h>
30 #include <TGeometry.h>
31 #include <TObjArray.h>
32 #include <TClonesArray.h>
36 #include <TParticle.h>
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"
46 #include "AliESDPmdTrack.h"
47 #include "AliESDEvent.h"
50 ClassImp(AliPMDtracker)
52 AliPMDtracker::AliPMDtracker():
54 fRecpoints(new TClonesArray("AliPMDrecpoint1", 1000)),
55 fPMDcontin(new TObjArray()),
56 fPMDcontout(new TObjArray()),
57 fPMDutil(new AliPMDUtility()),
69 // Default Constructor
72 //--------------------------------------------------------------------//
73 AliPMDtracker:: AliPMDtracker(const AliPMDtracker & /* tracker */):
74 TObject(/* tracker */),
91 AliError("Copy constructor not allowed");
94 //--------------------------------------------------------------------//
95 AliPMDtracker& AliPMDtracker::operator=(const AliPMDtracker & /* tracker */)
97 // assignment operator
98 AliError("Assignment operator not allowed");
102 //--------------------------------------------------------------------//
103 AliPMDtracker::~AliPMDtracker()
108 fRecpoints->Delete();
114 fPMDcontin->Delete();
120 fPMDcontout->Delete();
126 //--------------------------------------------------------------------//
127 void AliPMDtracker::LoadClusters(TTree *treein)
129 // Load the Reconstructed tree
132 //--------------------------------------------------------------------//
133 void AliPMDtracker::Clusters2Tracks(AliESDEvent *event)
135 // Converts digits to recpoints after running clustering
136 // algorithm on CPV plane and PREshower plane
143 TBranch *branch = fTreeR->GetBranch("PMDRecpoint");
146 AliError("PMDRecpoint branch not found");
149 branch->SetAddress(&fRecpoints);
151 Int_t nmodules = (Int_t) branch->GetEntries();
153 AliDebug(1,Form("Number of modules filled in treeR = %d",nmodules));
154 for (Int_t imodule = 0; imodule < nmodules; imodule++)
156 branch->GetEntry(imodule);
157 Int_t nentries = fRecpoints->GetLast();
158 AliDebug(2,Form("Number of clusters per modules filled in treeR = %d"
160 for(Int_t ient = 0; ient < nentries+1; ient++)
162 fPMDrecpoint = (AliPMDrecpoint1*)fRecpoints->UncheckedAt(ient);
163 idet = fPMDrecpoint->GetDetector();
164 ismn = fPMDrecpoint->GetSMNumber();
165 clusdata[0] = fPMDrecpoint->GetClusX();
166 clusdata[1] = fPMDrecpoint->GetClusY();
167 clusdata[2] = fPMDrecpoint->GetClusADC();
168 clusdata[3] = fPMDrecpoint->GetClusCells();
169 clusdata[4] = fPMDrecpoint->GetClusSigmaX();
170 clusdata[5] = fPMDrecpoint->GetClusSigmaY();
172 fPMDclin = new AliPMDrecpoint1(idet,ismn,clusdata);
173 fPMDcontin->Add(fPMDclin);
177 AliPMDDiscriminator *pmddiscriminator = new AliPMDEmpDiscriminator();
178 pmddiscriminator->Discrimination(fPMDcontin,fPMDcontout);
180 const Float_t kzpos = 361.5; // middle of the PMD
184 Float_t adc, ncell, rad;
185 Float_t xglobal = 0., yglobal = 0., zglobal = 0;
189 Int_t nentries2 = fPMDcontout->GetEntries();
190 AliDebug(1,Form("Number of clusters coming after discrimination = %d"
192 for (Int_t ient1 = 0; ient1 < nentries2; ient1++)
194 fPMDclout = (AliPMDclupid*)fPMDcontout->UncheckedAt(ient1);
196 det = fPMDclout->GetDetector();
197 smn = fPMDclout->GetSMN();
198 xpos = fPMDclout->GetClusX();
199 ypos = fPMDclout->GetClusY();
200 adc = fPMDclout->GetClusADC();
201 ncell = fPMDclout->GetClusCells();
202 rad = fPMDclout->GetClusRadius();
203 pid = fPMDclout->GetClusPID();
206 /**********************************************************************
207 * det : Detector, 0: PRE & 1:CPV *
208 * smn : Serial Module Number 0 to 23 for each plane *
209 * xpos : x-position of the cluster *
210 * ypos : y-position of the cluster *
211 * THESE xpos & ypos are not the true xpos and ypos *
212 * for some of the unit modules. They are rotated. *
213 * adc : ADC contained in the cluster *
214 * ncell : Number of cells contained in the cluster *
215 * rad : radius of the cluster (1d fit) *
216 **********************************************************************/
219 fPMDutil->RectGeomCellPos(smn,xpos,ypos,xglobal,yglobal);
223 zglobal = kzpos + 1.6; // PREshower plane
227 zglobal = kzpos - 1.7; // CPV plane
232 AliESDPmdTrack *esdpmdtr = new AliESDPmdTrack();
234 esdpmdtr->SetDetector(det);
235 esdpmdtr->SetClusterX(xglobal);
236 esdpmdtr->SetClusterY(yglobal);
237 esdpmdtr->SetClusterZ(zglobal);
238 esdpmdtr->SetClusterADC(adc);
239 esdpmdtr->SetClusterCells(ncell);
240 esdpmdtr->SetClusterPID(pid);
242 event->AddPmdTrack(esdpmdtr);
245 //--------------------------------------------------------------------//
246 void AliPMDtracker::SetVertex(Double_t vtx[3], Double_t evtx[3])
255 //--------------------------------------------------------------------//
256 void AliPMDtracker::ResetClusters()
258 if (fRecpoints) fRecpoints->Clear();
260 //--------------------------------------------------------------------//