]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TRD/qaRec/AliTRDcalibration.cxx
extend user interface
[u/mrichter/AliRoot.git] / TRD / qaRec / AliTRDcalibration.cxx
1
2 /**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 *                                                                        *
5 * Author: The ALICE Off-line Project.                                    *
6 * Contributors are mentioned in the code where appropriate.              *
7 *                                                                        *
8 * Permission to use, copy, modify and distribute this software and its   *
9 * documentation strictly for non-commercial purposes is hereby granted   *
10 * without fee, provided that the above copyright notice appears in all   *
11 * copies and that both the copyright notice and this permission notice   *
12 * appear in the supporting documentation. The authors make no claims     *
13 * about the suitability of this software for any purpose. It is          *
14 * provided "as is" without express or implied warranty.                  *
15 **************************************************************************/
16
17 // macro for very simple analysis
18
19
20
21 #include "Riostream.h"
22 #include "TChain.h"
23 #include "TTree.h"
24 #include "TProfile2D.h"
25 #include "TH2I.h"
26 #include "TH1F.h"
27 #include "TList.h"
28 #include "TMath.h"
29 #include "TCanvas.h"
30 #include "TObject.h"
31 #include "TFile.h"
32 #include "TObjArray.h"
33
34 #include "AliTRDrecoTask.h"
35 #include "AliAnalysisManager.h"
36
37 #include "AliESDInputHandler.h"
38 #include "AliTRDtrackV1.h"
39 #include "AliTRDseedV1.h"
40 #include "AliTRDcluster.h"
41 #include "AliTRDtrackInfo/AliTRDtrackInfo.h"
42 #include "AliTRDcalibDB.h"
43
44 #include "AliTRDCalibraFillHisto.h"
45
46 #include "AliLog.h"
47
48 #include "AliTRDcalibration.h"
49
50
51 ClassImp(AliTRDcalibration)
52
53 //________________________________________________________________________
54 AliTRDcalibration::AliTRDcalibration() 
55   :AliTRDrecoTask("Calibration", "Calibration on tracks")
56   ,fTrackInfo(0)
57   ,ftrdTrack(0)
58   ,fcl(0)
59   ,fTRDCalibraFillHisto(0)
60   ,fNbTRDTrackUsed(0)
61   ,fNbTimeBin(0x0)
62   ,fNbClusters(0)
63   ,fPHSum(0)
64   ,fCHSum(0)
65   ,flow(0)
66   ,fhigh(30)
67   ,fNbTimeBins(30)
68   ,ffillZero(kFALSE)
69 {
70   // Constructor
71 }  
72
73 //________________________________________________________________________
74 void AliTRDcalibration::CreateOutputObjects() 
75 {
76   OpenFile(0, "RECREATE");
77
78   // Number of time bins
79   AliTRDcalibDB *cal = AliTRDcalibDB::Instance();
80   fNbTimeBins = cal->GetNumberOfTimeBins();
81   
82   // instance calibration
83   fTRDCalibraFillHisto = AliTRDCalibraFillHisto::Instance();
84   fTRDCalibraFillHisto->SetHisto2d(); // choose to use histograms
85   fTRDCalibraFillHisto->SetCH2dOn();  // choose to calibrate the gain
86   fTRDCalibraFillHisto->SetPH2dOn();  // choose to calibrate the drift velocity
87   fTRDCalibraFillHisto->SetPRF2dOn(); // choose to look at the PRF
88   fTRDCalibraFillHisto->Init2Dhistos(); // initialise the histos
89   fTRDCalibraFillHisto->SetNumberClusters(flow); // At least 11 clusters
90   fTRDCalibraFillHisto->SetNumberClustersf(fhigh); // At least 11 clusters
91   fTRDCalibraFillHisto->SetFillWithZero(ffillZero); // Fill zeros
92
93   fTRDCalibraFillHisto->SetDebugLevel(fDebugLevel); //debug stuff
94
95   fContainer = new TObjArray();
96   fContainer->Add(fTRDCalibraFillHisto->GetCH2d()); //TH2I
97   fContainer->Add(fTRDCalibraFillHisto->GetPH2d()); //TProfile2D
98   fContainer->Add(fTRDCalibraFillHisto->GetPRF2d()); //TProfile2D
99
100   //printf("create output objects 1\n");
101
102   fNbTRDTrackUsed = new TH1F("TRDTrackUsed","TRDTrackUsed",50,0,50);
103   fNbTRDTrackUsed->Sumw2();
104   //
105   fNbTimeBin = new TH1F("TimeBin","TimeBin",35,0,35);
106   fNbTimeBin->Sumw2();
107   //
108   fNbClusters = new TH1F("NbClusters","",35,0,35);
109   fNbClusters->Sumw2();
110   //
111   fPHSum = new TProfile2D("PH2dSum","Nz0Nrphi0"
112       ,30,-0.05,(Double_t)fNbTimeBins/10.0-0.05
113       ,540,0,540);
114   fPHSum->SetYTitle("Det/pad groups");
115   fPHSum->SetXTitle("time [#mus]");
116   fPHSum->SetZTitle("<PH> [a.u.]");
117   fPHSum->SetStats(0);
118   //
119   fCHSum = new TH2I("CH2dSum","Nz0Nrphi0",100,0,300,540,0,540);
120   fCHSum->SetYTitle("Det/pad groups");
121   fCHSum->SetXTitle("charge deposit [a.u]");
122   fCHSum->SetZTitle("counts");
123   fCHSum->SetStats(0);
124   fCHSum->Sumw2();
125
126   fContainer->Add(fNbTRDTrackUsed);
127   fContainer->Add(fNbTimeBin);
128   fContainer->Add(fNbClusters);
129   fContainer->Add(fPHSum);
130   fContainer->Add(fCHSum);
131   //printf("create output objects 2\n");
132   
133 }
134
135 //________________________________________________________________________
136 void AliTRDcalibration::Exec(Option_t *) 
137 {
138   
139   Int_t nbTrdTracksUsed = 0;
140   
141   for(Int_t itrk=0; itrk < fTracks->GetEntriesFast(); itrk++){
142     
143     fTrackInfo = (AliTRDtrackInfo*)fTracks->UncheckedAt(itrk);
144     ftrdTrack = fTrackInfo->GetTrack();
145     if(!ftrdTrack) continue;
146     nbTrdTracksUsed++;
147     fTRDCalibraFillHisto->UpdateHistogramsV1(ftrdTrack);
148     
149     const AliTRDseedV1 *tracklet = 0x0;
150     for(Int_t itr = 0; itr < AliTRDgeometry::kNlayer; itr++){
151       if(!(tracklet = ftrdTrack->GetTracklet(itr))) continue;
152       if(!tracklet->IsOK()) continue;
153       Int_t nbclusters = 0;
154       // For PH
155       Double_t phtb[AliTRDseedV1::kNtb];
156       memset(phtb, 0, AliTRDseedV1::kNtb*sizeof(Double_t));
157
158       // For CH
159       Double_t sum = 0.0;
160       // normalisation
161       Float_t normalisation = 6.67;
162       Int_t detector = 0;
163       for(int ic=0; ic<AliTRDseedV1::kNTimeBins; ic++){
164         if(!(fcl = tracklet->GetClusters(ic))) continue;
165         nbclusters++;
166         Int_t time = fcl->GetPadTime();
167         Float_t ch =  tracklet->GetdQdl(ic);
168         detector = fcl->GetDetector();    
169         if((time>-1) && (time<fNbTimeBins)) phtb[time]=ch/normalisation;
170         sum += ch/normalisation;
171         //printf("time %d\n",time);
172         fNbTimeBin->Fill(time);
173       }
174       fNbClusters->Fill(nbclusters);
175       if((nbclusters > flow) && (nbclusters < fhigh)){
176         fCHSum->Fill(sum/20.0,0.0);
177         for(int ic=0; ic<fNbTimeBins; ic++){
178           if(ffillZero) fPHSum->Fill((Double_t)ic/10.0,0.0,(Double_t)phtb[ic]);
179           else {
180             if(phtb[ic] > 0.0) fPHSum->Fill((Double_t)ic/10.0,0.0,(Double_t)phtb[ic]);
181           }
182         }
183       }
184     }
185   }
186   
187   //Fill Histos
188   fNbTRDTrackUsed->Fill(nbTrdTracksUsed);
189   
190   
191   // Post output data
192   PostData(0, fContainer);
193   
194 }      
195 //________________________________________________________________________
196 void AliTRDcalibration::Terminate(Option_t *) 
197 {
198   // Draw result to the screen
199   // Called once at the end of the query
200
201   //printf("terminate\n");
202
203   if(fTRDCalibraFillHisto) fTRDCalibraFillHisto->DestroyDebugStreamer();
204
205   /*
206   fContainer = dynamic_cast<TList*> (GetOutputData(0));
207   if (!fContainer) {
208     Printf("ERROR: fList not available");
209     return;
210   }
211   */
212 }