]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/AliTRDtrackingResolution.cxx
Add resolution class
[u/mrichter/AliRoot.git] / TRD / qaRec / AliTRDtrackingResolution.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 /* $Id: AliTRDtrackingResolution.cxx 27496 2008-07-22 08:35:45Z cblume $ */
17
18 ////////////////////////////////////////////////////////////////////////////
19 //                                                                        //
20 //  Reconstruction QA                                                     //
21 //                                                                        //
22 //  Authors:                                                              //
23 //    Markus Fasel <M.Fasel@gsi.de>                                       //
24 //                                                                        //
25 ////////////////////////////////////////////////////////////////////////////
26
27 #include <TObjArray.h>
28 #include <TList.h>
29 #include <TH1F.h>
30 #include <TProfile.h>
31 #include <TMath.h>
32
33 #include "AliAnalysisManager.h"
34 #include "AliTRDseedV1.h"
35 #include "AliTrackReference.h"
36 #include "TTreeStream.h"
37
38 #include "AliTRDtrackInfo/AliTRDtrackInfo.h"
39 #include "AliTRDtrackingResolution.h"
40
41 #include <cstring>
42
43 ClassImp(AliTRDtrackingResolution)
44
45 //________________________________________________________
46 AliTRDtrackingResolution::AliTRDtrackingResolution(const char * name):
47   AliAnalysisTask(name, ""),
48   fTrackObjects(0x0),
49   fOutputHistograms(0x0),
50   fYres(0x0),
51 /*      fZres(0),
52   fYresAngle(0),
53   fPhiRes(0x0),
54   fPhiResAngle(0x0),*/
55   fDebugLevel(0),
56   fDebugStream(0x0)
57 {
58 //      memset(fYresLayer, 0, sizeof(TH1F *) * kNLayers);
59 //      memset(fZresLayer, 0, sizeof(TH1F *) * kNLayers);
60 //      memset(fYresLayerAngle, 0, sizeof(TProfile *) * kNLayers);
61 //      memset(fPhiResLayer, 0, sizeof(TH1F *) * kNLayers);
62 //      memset(fPhiResLayerAngle, 0, sizeof(TProfile *) * kNLayers);
63
64   DefineInput(0, TObjArray::Class());
65   DefineOutput(0, TList::Class());
66 }
67
68 //________________________________________________________
69 void AliTRDtrackingResolution::ConnectInputData(Option_t *){
70   fTrackObjects = dynamic_cast<TObjArray *>(GetInputData(0));
71 }
72
73 //________________________________________________________
74 void AliTRDtrackingResolution::CreateOutputObjects()
75 {
76   // spatial resolution
77   printf("Creating Histograms\n");
78   OpenFile(0, "RECREATE");
79   fOutputHistograms = new TList();
80   fYres = new TH1F("fYres", "y-Resolution", 100, -1.5, 1.5);
81   fOutputHistograms->Add(fYres);
82 //      fZres = new TH1F("fZres", "z-Resolution", 100, -1.5, 1.5);
83 //      fOutputHistograms->Add(fZres);
84 // 
85 //      fYresAngle = new TProfile("fYresAngle", "y-Resolution - Angluar dependence", 80, -40, 40);
86 //      fOutputHistograms->Add(fYresAngle);
87 // 
88 //      // angular resolution
89 //      fPhiRes = new TH1F("fPhiRes", "phi-resolution", 20, -10, 10);
90 //      fOutputHistograms->Add(fPhiRes);
91 //      
92 //      fPhiResAngle = new TProfile("fPhiResAngle", "phi-resolution - Angular dependence", 80, -40, 40);
93 //      fOutputHistograms->Add(fPhiResAngle);
94
95 /*      for(Int_t iplane = 0; iplane < kNLayers; iplane++){
96     // spatial resolution
97     fYresLayer[iplane] = new TH1F(Form("fYresLayer%d", iplane), Form("y-Resolution in Layer %d", iplane), 100, -1.5, 1.5);
98     fOutputHistograms->Add(fYresLayer[iplane]);
99
100     fZresLayer[iplane] = new TH1F(Form("fZresLayer%d", iplane), Form("z-Resolution in Layer %d", iplane), 100, -1.5, 1.5);
101     fOutputHistograms->Add(fZresLayer[iplane]);
102
103     fYresLayerAngle[iplane] = new TProfile(Form("fYresLayerAngle%d", iplane), Form("y-Resolution in Layer %d - Angluar dependence", iplane), 80, -40, 40);
104     fOutputHistograms->Add(fYresLayerAngle[iplane]);
105
106     // angular resolution
107     fPhiResLayer[iplane] = new TH1F(Form("fPhiResLayer%d", iplane), Form("phi-resolution in Layer %d", iplane), 20, -10, 10);
108     fOutputHistograms->Add(fPhiResLayer[iplane]);
109
110     fPhiResLayerAngle[iplane] = new TProfile(Form("fPhiResAngle%d", iplane), Form("phi-resolution in Layer %d - Angular dependence", iplane), 80, -40, 40);
111     fOutputHistograms->Add(fPhiResLayerAngle[iplane]);
112   }*/
113 }
114
115 //________________________________________________________
116 void AliTRDtrackingResolution::Exec(Option_t *){
117   // spatial Resolution: res = pos_{Tracklet}(x = x_{Anode wire}) - pos_{TrackRef}(x = x_{Anode wire})
118   // angular Resolution: res = Tracklet angle - TrackRef Angle
119   Int_t nTrackInfos = fTrackObjects->GetEntriesFast();
120   if(fDebugLevel>=2) printf("Number of Histograms: %d\n", fOutputHistograms->GetEntries());
121   AliTRDtrackInfo *fInfo = 0x0;
122   if(fDebugLevel>=2) printf("Number of TrackInfos: %d\n", nTrackInfos);
123   for(Int_t iTI = 0; iTI < nTrackInfos; iTI++){
124     if(fDebugLevel>=2) printf("Doing Object %d\n", iTI);
125     fInfo = dynamic_cast<AliTRDtrackInfo *>(fTrackObjects->UncheckedAt(iTI));
126     // check if ESD and MC-Information are available
127     if(!fInfo || !fInfo->GetTRDtrack() || fInfo->GetNTrackRefs() < 2) continue; 
128     AliTRDseedV1*fTracklet = 0;
129     AliTrackReference *fTrackRefs[2];
130     for(Int_t iplane = 0; iplane < kNLayers; iplane++){
131       if(fDebugLevel>=2) printf("plane %d\n", iplane);
132       fTracklet = fInfo->GetTracklet(iplane);
133       if(!fTracklet) continue;
134       // check for 2 track ref where the radial position has a distance less than 3.7mm
135       if(fDebugLevel>=2) printf("Find Track References for x = %f\n", fTracklet->GetX0());
136       if(fDebugLevel>=2) printf("Number of Clusters: %d\n", fTracklet->GetN());
137       Int_t nFound = 0;
138       memset(fTrackRefs, 0, sizeof(AliTrackReference*) * 2);
139       AliTrackReference *tempTrackRef = 0;
140       for(Int_t itr = 0; itr < fInfo->GetNTrackRefs(); itr++){
141         if(fDebugLevel>=2) printf("nFound = %d\n", nFound);
142         if(nFound >= 2) break;
143         tempTrackRef = fInfo->GetTrackRef(itr);
144         if(!tempTrackRef) continue;
145         if(fDebugLevel>=2) printf("TrackRef %d: x = %f\n", itr, tempTrackRef->LocalX());
146         if(fTracklet->GetX0() - tempTrackRef->LocalX() > 3.7) continue;
147         if(tempTrackRef->LocalX() - fTracklet->GetX0() > 3.7) break;
148         if(fDebugLevel>=2) printf("accepted\n");
149         if(nFound == 1)
150           if(fTrackRefs[0]->LocalX() >= tempTrackRef->LocalX()) continue;
151         fTrackRefs[nFound++] = tempTrackRef;
152       }
153       if(fDebugLevel>=2) printf("nFound = %d\n", nFound);
154       if(nFound < 2) continue;
155       // We found 2 track refs for the tracklet, get y and z at the anode wire by a linear approximation
156       Double_t dx = fTrackRefs[1]->LocalX() - fTrackRefs[0]->LocalX();
157       Double_t dydx = (fTrackRefs[1]->LocalY() - fTrackRefs[0]->LocalY()) / dx;
158       Double_t dzdx = (fTrackRefs[1]->Z() - fTrackRefs[0]->Z()) / dx;
159       Double_t dx0 = fTrackRefs[1]->LocalX() - fTracklet->GetX0();
160       Double_t ymc =  fTrackRefs[1]->LocalY() - dydx*dx0;
161       Double_t zmc =  fTrackRefs[1]->Z() - dzdx*dx0;
162       
163       Double_t dy = fTracklet->GetYfit(0) - ymc;
164       Double_t dz = fTracklet->GetZfit(0) - zmc;
165       //res_y *= 100; // in mm
166       Double_t momentum = fTrackRefs[0]->P();
167       
168       // Fill Histograms
169       if(fDebugLevel>=2) printf("dy = %f\n", dy);
170       fYres->Fill(dy);
171 //                      fZres->Fill(res_z);
172 /*                      fYresLayer[iplane]->Fill(res_y);
173       fZresLayer[iplane]->Fill(res_z);*/
174
175 //       Double_t phi     = fTrackRefs[0]->Phi();
176 //       Double_t theta   = fTrackRefs[0]->Theta();
177       Double_t phi   = TMath::ATan(dydx);
178       Double_t theta = TMath::ATan(dzdx);
179       
180 //                      fYresAngle->Fill(phi, res_y);
181 //                      fYresLayerAngle[iplane]->Fill(phi, res_y);
182       
183       Double_t dphi   = TMath::ATan(fTracklet->GetZfit(1)) - phi;
184       
185 //                      fPhiRes->Fill(dphi);
186 //                      fPhiResLayer[iplane]->Fill(dphi);
187 //                      fPhiResAngle->Fill(phi, dphi);
188 //                      fPhiResLayerAngle[iplane]->Fill(phi, dphi);
189       
190       // Fill Debug Tree
191       if(fDebugLevel>=1){
192         (*fDebugStream) << "Resolution"
193           << "plane="           << iplane
194           << "p="       << momentum
195           << "dx="      << dx
196           << "dy="                << dy
197           << "dz="                << dz
198           << "phi="                     << phi
199           << "theta="           << theta
200           << "dphi="            << dphi
201           << "\n";
202       }
203     }
204   }
205   PostData(0, fOutputHistograms);
206 }
207
208 //________________________________________________________
209 void AliTRDtrackingResolution::Terminate(Option_t *){
210   //printf("Tracking Resolution: Terminate\n");
211   if(fDebugStream) delete fDebugStream;
212 }
213
214 //________________________________________________________
215 void AliTRDtrackingResolution::SetDebugLevel(Int_t level){
216   fDebugLevel = level;
217   if(!fDebugLevel) return;
218   if(fDebugStream) return;
219   fDebugStream = new TTreeSRedirector("TRD.Resolution.root");
220 }