]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TRD/AliTRDmultiplicity.cxx
second interation on implementing variable detector segmentation
[u/mrichter/AliRoot.git] / PWG1 / TRD / AliTRDmultiplicity.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 //                                                                        //
19 //  TRD multiplicity                                                      //
20 //                                                                        //
21 // Task to select true single tracks
22 // 
23 // Program flow:
24 // Part A - task AliTRDmultiplicity (within the framework qaRec):
25 //   For TRD standalone or TRD&TOF tracks I make a FitRiemanTilt
26 //   The resulting parameterization is dumped into a DebugStreamer written to a file
27 // 
28 // Part B – $ALICE_ROOT/TRD/qaRec/macros/TrackletsinTRD.C[h] (analysis macro):
29 //   The clusters are read in
30 //   The above produced file is read in
31 //   I make a fit through the parameterization of the  FitRiemanTilt -> in order to get a straight line (representing the track)
32 //   the distance of each cluster with respect to the straight line is calculated
33 //   If there is no cluster at a distance of 0.6  to 2 with respect to the track, the track is defined as a good track, i.e. clean single track
34 //                                                                        //
35 //  Authors: Yvonne C. Pachmayer <pachmay@physi.uni-heidelberg.de>        //
36 //                                                                        //
37 ////////////////////////////////////////////////////////////////////////////
38
39 #include <TObjArray.h>
40 #include <TH1F.h>
41 #include <TTreeStream.h>
42
43 #include "AliESDtrack.h"
44 #include "AliExternalTrackParam.h"
45 #include "AliTRDseedV1.h"
46 #include "AliTRDtrackV1.h"
47 #include "AliTRDtrackerV1.h"
48 #include "AliTrackPointArray.h"
49
50 #include "AliTRDmultiplicity.h"
51 #include "AliTRDgeometry.h"
52 #include "info/AliTRDtrackInfo.h"
53
54 ClassImp(AliTRDmultiplicity)
55
56 //____________________________________________________________________
57 AliTRDmultiplicity::AliTRDmultiplicity()
58     :AliTRDrecoTask()
59   ,fEventCounter(0)
60 {
61   //
62   // Default constructor
63     //
64 }
65
66 AliTRDmultiplicity::AliTRDmultiplicity(char* name)
67     :AliTRDrecoTask(name, "Barrel Tracking Multiplicity")
68   ,fEventCounter(0)
69 {
70   //
71   // Default constructor
72     //
73 }
74
75 //____________________________________________________________________
76 AliTRDmultiplicity::~AliTRDmultiplicity()
77 {
78 }
79
80 //____________________________________________________________________
81 void  AliTRDmultiplicity::UserCreateOutputObjects()
82 {
83   //
84   // Create output objects
85   //
86
87   OpenFile(1, "RECREATE");
88
89   TH1 *h = 0x0;
90   fContainer = new TObjArray();
91   for(Int_t is=0; is<AliTRDgeometry::kNsector; is++){
92   fContainer->Add(h = new TH1F(Form("h_sector_%i", is), " ", 100,-10,10));
93   }
94
95
96 //____________________________________________________________________
97 void AliTRDmultiplicity::UserExec(Option_t *)
98 {
99   //
100   // Do it
101   //
102
103
104
105   ULong_t status;
106   AliTRDtrackInfo     *track = 0x0;
107   AliExternalTrackParam *esd = 0x0;
108   AliTRDtrackV1 *trackTRD = 0x0;
109   Double_t x_anode[AliTRDgeometry::kNlayer] = {300.2, 312.8, 325.4, 338.0, 350.6, 363.2}; // Take the default X0
110   Int_t standAlone=1;
111   fEventCounter++;
112
113
114   for(Int_t itrk=0; itrk<fTracks->GetEntriesFast(); itrk++){
115     track = (AliTRDtrackInfo*)fTracks->UncheckedAt(itrk);
116     if(!track->HasESDtrack()) continue;
117     status = track->GetStatus();
118
119     AliTrackPoint points[6];
120     Float_t xyz[3];
121     memset(xyz, 0, sizeof(Float_t) * 3);
122     Float_t trackletX[6];
123     Float_t trackletY[6];
124     Float_t trackletZ[6];
125     for(Int_t a=0;a<6;a++)
126     {
127         xyz[0] = x_anode[a];
128         points[a].SetXYZ(xyz);
129         trackletX[a]=-1000;
130         trackletY[a]=-1000;
131         trackletZ[a]=-1000;
132     }
133     Int_t detTracklet=600;
134
135
136
137     // TRD standalone
138     if(((status&AliESDtrack::kTRDout)>0) && !((status&AliESDtrack::kTRDin)>0))
139     {
140         trackTRD = track->GetTrack();
141
142         if(trackTRD)
143         {
144             AliTRDseedV1 *tracklet = 0x0;
145             Int_t counterTracklet=0;
146             for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++)
147             {
148                 tracklet = trackTRD->GetTracklet(itl);
149                 if(tracklet)
150                 {
151                     if(tracklet->IsOK())
152                     {
153                         counterTracklet++;
154                         detTracklet=tracklet->GetDetector();
155                         //printf("%i %f %f \n",itl,tracklet->GetYfit(0),tracklet->GetZfit(0));
156                         trackletX[itl]=tracklet->GetX0();
157                         trackletY[itl]=tracklet->GetYfit(0);
158                         trackletZ[itl]=tracklet->GetZfit(0);
159                     }
160                 }
161             }
162             // this cut is needed otherwise the order of tracklets in the fit function is not correct
163             if(counterTracklet==AliTRDgeometry::kNlayer) AliTRDtrackerV1::FitRiemanTilt(const_cast<AliTRDtrackV1 *>(trackTRD), 0x0, kTRUE, counterTracklet, points);
164             else continue;
165
166             if(DebugLevel()>=1){
167               for(Int_t il = 0; il < AliTRDgeometry::kNlayer; il++){
168                 //printf("---------------------------------------- %i %i %f %f %f \n",counterTracklet,il,points[il].GetX(),points[il].GetY(),points[il].GetZ());
169                 Double_t pointX=points[il].GetX();
170                 Double_t pointY=points[il].GetY();
171                 Double_t pointZ=points[il].GetZ();
172
173
174                 (*DebugStream())   << "TrackletsinTRD"
175                     << "standalone=" << standAlone
176                     << "eventcounter=" << fEventCounter
177                     << "layer="     << il
178                     << "dettracklet=" << detTracklet
179                     << "xtracklet=" << trackletX[il]
180                     << "xtrack="    << pointX
181                     << "ytracklet=" << trackletY[il]
182                     << "ytrack="    << pointY
183                     << "ztracklet=" << trackletZ[il]
184                     << "ztrack="    << pointZ
185                     << "num_tracklets=" << counterTracklet
186                     << "\n";
187
188               }
189             }
190         }
191     } // end TRD standalone
192
193
194
195     // TPC cluster selection for cosmic data
196     if((track->GetTPCncls())<40) continue;
197     // missing TPC propagation
198     if(!(status&AliESDtrack::kTPCout)) continue;
199     standAlone=0;
200
201     esd = track->GetESDinfo()->GetOuterParam();
202     // calculate sector - does not matter if track ends in TPC or TRD - sector number independent
203     Double_t alpha;
204     alpha=-100;
205     Int_t fSector=100;
206     if(esd){
207         alpha=esd->GetAlpha();
208         fSector = (Int_t)((alpha+ TMath::Pi())/(20*TMath::Pi()/180));
209         if((fSector>-1) && (fSector<18))
210         {
211         }
212         else
213         {
214             continue;
215         }
216     }
217
218
219
220     // only these sectors have a TRD detector at the moment
221    if(fSector==0||fSector==8||fSector==9||fSector==17)
222     {
223         trackTRD = track->GetTrack();
224         if(!trackTRD) continue;
225         AliTRDseedV1 *tracklet = 0x0;
226         Int_t counterTracklet=0;
227
228
229         for(Int_t itl = 0; itl < AliTRDgeometry::kNlayer; itl++){
230             tracklet = trackTRD->GetTracklet(itl);
231             if(!tracklet || !tracklet->IsOK()) continue;
232             counterTracklet++;
233             detTracklet=tracklet->GetDetector();
234             trackletX[itl]=tracklet->GetX0();
235             trackletY[itl]=tracklet->GetYfit(0);
236             trackletZ[itl]=tracklet->GetZfit(0);
237
238             /*
239             AliTRDcluster *c = 0x0;
240             for(Int_t itime = 0; itime < AliTRDseedV1::kNtb; itime++){
241                 c = tracklet->GetClusters(itime);
242                 if(!c) continue;
243 //                Float_t cluster_x = c->GetX();
244 //                Bool_t isprop;
245 //                isprop=kFALSE;
246 //                isprop=AliTracker::PropagateTrackTo(esd,(Double_t)cluster_x,0.105,5,kFALSE);
247 //                if(isprop)
248 //                {
249 //                    Int_t detector = c->GetDetector();
250 //                    ((TH1F*)fContainer->At(fSector))->Fill((c->GetY())-(esd->GetY()));
251 //                    if(c->GetY()-esd->GetY()>);
252 //                    printf("diff: %f\n",c->GetY()-(esd->GetY()));
253 //                }
254             }
255             */
256
257         } // loop over tracklets ends
258
259         if(counterTracklet==AliTRDgeometry::kNlayer)
260         {
261             AliTRDtrackerV1::FitRiemanTilt(const_cast<AliTRDtrackV1 *>(trackTRD), 0x0, kTRUE, counterTracklet, points);
262         }
263         else continue;
264
265        
266
267         if(DebugLevel()>=1){
268           for(Int_t il = 0; il < AliTRDgeometry::kNlayer; il++){
269             //printf("---------------------------------------- %i %i %f %f %f \n",counterTracklet,il,points[il].GetX(),points[il].GetY(),points[il].GetZ());
270             Double_t pointX=points[il].GetX();
271             Double_t pointY=points[il].GetY();
272             Double_t pointZ=points[il].GetZ();
273             
274             (*DebugStream())   << "TrackletsinTRD"
275                 << "standalone=" << standAlone
276                 << "eventcounter=" << fEventCounter
277                 << "layer="     << il
278                 << "dettracklet=" << detTracklet
279                 << "xtracklet=" << trackletX[il]
280                 << "xtrack="    << pointX
281                 << "ytracklet=" << trackletY[il]
282                 << "ytrack="    << pointY
283                 << "ztracklet=" << trackletZ[il]
284                 << "ztrack="    << pointZ
285                 << "num_tracklets=" << counterTracklet
286                 << "\n";
287           }
288         }
289
290     }
291   }
292
293   
294   PostData(1, fContainer);
295 }
296
297