]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/TPC/AliTPCPerformanceSummary.cxx
compilation worning removed
[u/mrichter/AliRoot.git] / PWG1 / TPC / AliTPCPerformanceSummary.cxx
1 //------------------------------------------------------------------------------
2 // Implementation of AliTPCPerformanceSummary class. 
3 // It has only static member functions to extract some TPC Performance
4 // parameters and produce trend graphs.
5 // The function MakeReport is to be called for every run. It reads AliPerformanceTPC
6 // and AliPerformanceDEdx objects from file and produces a
7 // rootfile with the results stored in a TTree.
8 // The function MakeReport needs a list of these rootfiles as input
9 // and writes the output (tree and histograms) to another rootfile.
10 //
11 // Author: M.Knichel 2010-05-21
12 //------------------------------------------------------------------------------
13
14 #include <fstream>
15
16 #include "TSystem.h"
17 #include "TMath.h"
18 #include "TVectorD.h"
19 #include "TList.h"
20 #include "TFile.h"
21 #include "TGrid.h"
22 #include "TF1.h"
23 #include "TH1.h"
24 #include "TH2.h"
25 #include "TProfile.h"
26 #include "THnSparse.h"
27 #include "TTree.h"
28 #include "TChain.h"
29 #include "TGraph.h"
30 #include "TPad.h"
31 #include "TCanvas.h"
32
33 #include "TTreeStream.h"
34 #include "AliPerformanceTPC.h"
35 #include "AliPerformanceDEdx.h"
36 #include "AliPerformanceDCA.h"
37
38 #include "AliTPCPerformanceSummary.h"
39
40 ClassImp(AliTPCPerformanceSummary)
41
42 //_____________________________________________________________________________
43 Int_t AliTPCPerformanceSummary::WriteToTTreeSRedirector(const AliPerformanceTPC* pTPC, const AliPerformanceDEdx* pTPCgain, TTreeSRedirector* pcstream, Int_t run)
44 {
45     // 
46     // Extracts performance parameters from pTPC and pTPCgain.
47     // Output is written to pcstream.
48     // The run number must be provided since it is not stored in 
49     // AliPerformanceTPC or AliPerformanceDEdx.
50     //
51     
52     if (!pcstream) return -1;
53     (*pcstream)<<"tpcQA"<<"run="<<run;
54     Int_t returncode = 0;
55     returncode += AnalyzeNCL(pTPC, pcstream);    
56     returncode += AnalyzeDrift(pTPC, pcstream);
57     returncode += AnalyzeDriftPos(pTPC, pcstream);
58     returncode += AnalyzeDriftNeg(pTPC, pcstream);
59     returncode += AnalyzeGain(pTPCgain, pcstream);
60     returncode += AnalyzeDCARPhi(pTPC, pcstream);
61     returncode += AnalyzeDCARPhiPos(pTPC, pcstream);
62     returncode += AnalyzeDCARPhiNeg(pTPC, pcstream);
63     returncode += AnalyzeEvent(pTPC, pcstream);
64     (*pcstream)<<"tpcQA"<<"\n";
65     return returncode;
66 }
67
68 //_____________________________________________________________________________
69 Int_t AliTPCPerformanceSummary::WriteToFile(const AliPerformanceTPC* pTPC, const AliPerformanceDEdx* pTPCgain, const Char_t* outfile, Int_t run)
70 {
71     //
72     // Extracts performance parameters from pTPC and pTPCgain.
73     // Output is written to a TTree saved in outfile.
74     // The run number must be provided since it is not stored in 
75     // AliPerformanceTPC or AliPerformanceDEdx.
76     //
77     // The function creates a TTreeSRedirector and calls the 
78     // function WriteToTTreeSRedirector.
79     //
80     
81     if (!outfile) return -1;
82     TTreeSRedirector* pcstream = 0;
83     pcstream = new TTreeSRedirector(outfile);
84     if (!pcstream) return -1;
85     Int_t returncode = WriteToTTreeSRedirector(pTPC, pTPCgain, pcstream, run);
86     if (pcstream) { delete pcstream; pcstream = 0; }
87     return returncode;
88     
89 }
90
91 //_____________________________________________________________________________
92 Int_t AliTPCPerformanceSummary::MakeReport(const Char_t* infile, const Char_t* outfile, Int_t run)
93 {
94     //
95     // Reads QA information (AliPerformanceTPC and AliPerformanceDEdx) from
96     // infile (this must be a rootfile) and writes the output to a TTree
97     // stored in outfile.
98     // The run number must be provided since it is not stored in 
99     // AliPerformanceTPC or AliPerformanceDEdx.
100     // 
101     // The input objects must be named "AliPerformanceTPC" and 
102     // "AliPerformanceDEdxTPCInner" and stored in a TList which name must
103     // be one of the following: "TPC", "TPCQA", "TPC_PerformanceQA"
104     // or "TPC_PerformanceQA/TPC" (with directory)
105     //
106     
107     if (!infile) return -1;
108     if (!outfile) return -1;
109     TFile *f =0;
110     f=TFile::Open(infile,"read");
111     if (!f) {
112         printf("File %s not available\n", infile);
113         return -1;
114     } 
115     TList* list = 0;
116     list = dynamic_cast<TList*>(f->Get("TPC")); 
117     if (!list) { list = dynamic_cast<TList*>(f->Get("TPCQA")); }
118     if (!list) { list = dynamic_cast<TList*>(f->Get("TPC_PerformanceQA")); }
119     if (!list) { list = dynamic_cast<TList*>(f->Get("TPC_PerformanceQA/TPCQA")); }
120     if (!list) {
121             printf("QA %s not available\n", infile);
122             return -1;
123     } 
124     AliPerformanceTPC* pTPC = 0;
125     AliPerformanceDEdx* pTPCgain = 0; 
126     if (list) {  pTPC = (AliPerformanceTPC*)(list->FindObject("AliPerformanceTPC")); }
127     if (list) {  pTPCgain = (AliPerformanceDEdx*)(list->FindObject("AliPerformanceDEdxTPCInner")); }
128     
129     Int_t returncode = 0;
130     returncode = WriteToFile(pTPC, pTPCgain, outfile, run);
131     if (f) { delete f; f=0; }
132     return returncode;
133 }
134
135 //_____________________________________________________________________________
136 Int_t AliTPCPerformanceSummary::ProduceTrends(const Char_t* infilelist, const Char_t* outfile)
137 {
138     //
139     // Produces trend graphs.
140     //
141     // Input: infilelist is a textfile with one rootfile per line.
142     // There should be one rootfile for each run, the rootfile must
143     // contain the output of the MakeReport function
144     // Output: the information for all runs is merged into a TTree
145     // that is saved in outfile along with the trend graphs.
146     // Trend graphs are stored as TCanvas objects to include axis labels etc.
147     //
148     
149     if (!infilelist) return -1;
150     if (!outfile) return -1;
151      
152     TChain* chain = new TChain("tpcQA");
153     ifstream in;
154     in.open(infilelist);
155
156     TString currentFile;    
157     while(in.good()) {
158         in >> currentFile;
159
160         if (!currentFile.Contains("root")) continue; // protection            
161         chain->Add(currentFile.Data());
162     }
163     in.close();
164     //TTree *tree = chain;
165     TTree *tree = chain->CopyTree("1");
166     if (chain) { delete chain; chain=0; }
167     //TGraph* graph = dynamic_cast<TGraph*>(tree->DrawClone("run:run"));
168     //TGraph *graph = (TGraph*)gPad->GetPrimitive("Graph");
169     
170     TFile* out = new TFile(outfile,"RECREATE");
171     out->cd();
172     Char_t* condition = "meanTPCncl>0";
173     SaveGraph(tree,"meanTPCnclF","run",condition);
174     SaveGraph(tree,"rmsTPCnclF","run",condition);
175     SaveGraph(tree,"meanTPCChi2","run",condition);
176     SaveGraph(tree,"rmsTPCChi2","run",condition);
177     SaveGraph(tree,"slopeATPCnclF","run",condition);
178     SaveGraph(tree,"slopeCTPCnclF","run",condition);
179     SaveGraph(tree,"slopeATPCnclFErr","run",condition);
180     SaveGraph(tree,"slopeCTPCnclFErr","run",condition);
181     SaveGraph(tree,"meanTPCncl","run",condition);
182     SaveGraph(tree,"rmsTPCncl","run",condition);
183     SaveGraph(tree,"slopeATPCncl","run",condition);
184     SaveGraph(tree,"slopeCTPCncl","run",condition);
185     SaveGraph(tree,"slopeATPCnclErr","run",condition);
186     SaveGraph(tree,"slopeCTPCnclErr","run",condition);
187     
188     SaveGraph(tree,"offsetdRA","run",condition);
189     SaveGraph(tree,"slopedRA","run",condition);
190     SaveGraph(tree,"offsetdRC","run",condition);
191     SaveGraph(tree,"slopedRC","run",condition);
192     SaveGraph(tree,"offsetdRAErr","run",condition);
193     SaveGraph(tree,"slopedRAErr","run",condition);    
194     SaveGraph(tree,"offsetdRCErr","run",condition);
195     SaveGraph(tree,"slopedRCErr","run",condition);
196     SaveGraph(tree,"offsetdRAchi2","run",condition);
197     SaveGraph(tree,"slopedRAchi2","run",condition);
198     SaveGraph(tree,"offsetdRCchi2","run",condition);    
199     SaveGraph(tree,"slopedRCchi2","run",condition);    
200     
201     SaveGraph(tree,"offsetdRAPos","run",condition);
202       
203     SaveGraph(tree,"slopedRAPos","run",condition);
204     SaveGraph(tree,"offsetdRCPos","run",condition);
205     SaveGraph(tree,"slopedRCPos","run",condition);
206     SaveGraph(tree,"offsetdRAErrPos","run",condition);
207     SaveGraph(tree,"slopedRAErrPos","run",condition);
208     SaveGraph(tree,"offsetdRCErrPos","run",condition); 
209     SaveGraph(tree,"slopedRCErrPos","run",condition);
210     SaveGraph(tree,"offsetdRAchi2Pos","run",condition);
211     SaveGraph(tree,"slopedRAchi2Pos","run",condition);
212     SaveGraph(tree,"offsetdRCchi2Pos","run",condition);
213     SaveGraph(tree,"slopedRCchi2Pos","run",condition);
214         
215     SaveGraph(tree,"offsetdRANeg","run",condition);
216     SaveGraph(tree,"slopedRANeg","run",condition);
217     SaveGraph(tree,"offsetdRCNeg","run",condition);
218     SaveGraph(tree,"slopedRCNeg","run",condition);
219     SaveGraph(tree,"offsetdRAErrNeg","run",condition);
220     SaveGraph(tree,"slopedRAErrNeg","run",condition);
221     SaveGraph(tree,"offsetdRCErrNeg","run",condition);
222     SaveGraph(tree,"slopedRCErrNeg","run",condition);
223     SaveGraph(tree,"offsetdRAchi2Neg","run",condition);
224     SaveGraph(tree,"slopedRAchi2Neg","run",condition);
225     SaveGraph(tree,"offsetdRCchi2Neg","run",condition);
226     SaveGraph(tree,"slopedRCchi2Neg","run",condition);
227         
228     SaveGraph(tree,"offsetdZAPos","run",condition);
229     SaveGraph(tree,"slopedZAPos","run",condition);
230     SaveGraph(tree,"offsetdZCPos","run",condition);
231     SaveGraph(tree,"slopedZCPos","run",condition);
232     SaveGraph(tree,"offsetdZAErrPos","run",condition);
233     SaveGraph(tree,"slopedZAErrPos","run",condition);
234     SaveGraph(tree,"offsetdZCErrPos","run",condition);
235     SaveGraph(tree,"slopedZCErrPos","run",condition);
236     SaveGraph(tree,"offsetdZAchi2Pos","run",condition);
237     SaveGraph(tree,"slopedZAchi2Pos","run",condition);
238     SaveGraph(tree,"offsetdZCchi2Pos","run",condition);
239     SaveGraph(tree,"slopedZCchi2Pos","run",condition);
240     
241     SaveGraph(tree,"offsetdZANeg","run",condition);
242     SaveGraph(tree,"slopedZANeg","run",condition);
243     SaveGraph(tree,"offsetdZCNeg","run",condition);
244     SaveGraph(tree,"slopedZCNeg","run",condition);
245     SaveGraph(tree,"offsetdZAErrNeg","run",condition);
246     SaveGraph(tree,"slopedZAErrNeg","run",condition);
247     SaveGraph(tree,"offsetdZCErrNeg","run",condition);
248     SaveGraph(tree,"slopedZCErrNeg","run",condition);
249     SaveGraph(tree,"offsetdZAchi2Neg","run",condition);
250     SaveGraph(tree,"slopedZAchi2Neg","run",condition);
251     SaveGraph(tree,"offsetdZCchi2Neg","run",condition);
252     SaveGraph(tree,"slopedZCchi2Neg","run",condition);    
253     
254     SaveGraph(tree,"offsetdZA","run",condition);
255     SaveGraph(tree,"slopedZA","run",condition);
256     SaveGraph(tree,"offsetdZC","run",condition);
257     SaveGraph(tree,"slopedZC","run",condition);
258     SaveGraph(tree,"offsetdZAErr","run",condition);
259     SaveGraph(tree,"slopedZAErr","run",condition);
260     SaveGraph(tree,"offsetdZCErr","run",condition);
261     SaveGraph(tree,"slopedZCErr","run",condition);
262     SaveGraph(tree,"offsetdZAchi2","run",condition);
263     SaveGraph(tree,"slopedZAchi2","run",condition);
264     SaveGraph(tree,"offsetdZCchi2","run",condition);
265     SaveGraph(tree,"slopedZCchi2","run",condition);        
266
267     SaveGraph(tree,"meanVertX","run",condition);
268     SaveGraph(tree,"rmsVertX","run",condition);
269     SaveGraph(tree,"meanVertY","run",condition);
270     SaveGraph(tree,"rmsVertY","run",condition);
271     SaveGraph(tree,"meanVertZ","run",condition);
272     SaveGraph(tree,"rmsVertZ","run",condition);
273     SaveGraph(tree,"vertStatus","run",condition);
274     SaveGraph(tree,"meanMult","run",condition);
275     SaveGraph(tree,"rmsMult","run",condition);
276     SaveGraph(tree,"meanMultPos","run",condition);
277     SaveGraph(tree,"rmsMultPos","run",condition);
278     SaveGraph(tree,"meanMultNeg","run",condition);
279     SaveGraph(tree,"rmsMultNeg","run",condition);
280     SaveGraph(tree,"vertAll","run",condition);
281     SaveGraph(tree,"vertOK","run",condition);
282
283     tree->Write();
284     
285     out->Close();   
286     if (tree) { delete tree; tree=0; }
287     if (out) { delete out; out=0; }
288 return 0;
289 }
290
291 //_____________________________________________________________________________
292 Int_t AliTPCPerformanceSummary::SaveGraph(TTree* tree, const Char_t* y, const Char_t* x, const Char_t* condition)
293 {    
294     //
295     // Creates a Graph and writes the canvas to the current directory
296     // called by ProduceTrends function.
297     //
298     
299     TString s(y);
300     s += ':';
301     s += x;
302     tree->Draw(s.Data(),condition,"goff");
303     TCanvas* c = new TCanvas(s.Data(),s.Data());
304     c->cd();
305     TPad* p = new TPad("pad0","pad0",0,0,1,1);
306     p->Draw();
307     p->cd();
308     if (tree->GetSelectedRows() > 0) {
309       TGraph* graph = new TGraph(tree->GetSelectedRows(), tree->GetV2(), tree->GetV1());
310       graph->Draw("AP");
311       graph->GetXaxis()->SetTitle(x);
312       graph->GetYaxis()->SetTitle(y);
313       c->Write(s.Data());
314       delete graph;
315     }
316     //graph->Write(s.Data());
317     delete c;
318     
319 return 0;
320 }
321
322 //_____________________________________________________________________________
323 Int_t AliTPCPerformanceSummary::AnalyzeDCARPhi(const AliPerformanceTPC* pTPC, TTreeSRedirector* pcstream)
324 {
325     //
326     // Analyse DCA R imperfections
327     //
328     
329     if (!pcstream) return 8;
330     if (!pTPC) return 8;
331         
332     // variables:
333     static Double_t offsetdRA=0;
334     static Double_t slopedRA=0;
335     static Double_t offsetdRC=0;
336     static Double_t slopedRC=0;
337     static Double_t offsetdRAErr=0;
338     static Double_t slopedRAErr=0;
339     static Double_t offsetdRCErr=0;
340     static Double_t slopedRCErr=0;
341     static Double_t offsetdRAchi2=0;
342     static Double_t slopedRAchi2=0;
343     static Double_t offsetdRCchi2=0;
344     static Double_t slopedRCchi2=0;
345
346     //AliPerformanceTPC* pTPC =  dynamic_cast<AliPerformanceTPC*>(pTPCObject);    
347     TH1* his1D=0;
348     TH2* his2D=0;
349     static TF1 *fpol1 = new TF1("fpol1","pol1");
350     TObjArray arrayFit;
351     his2D = pTPC->GetTPCTrackHisto()->Projection(3,5); 
352     his2D->FitSlicesY(0,0,-1,10,"QNR",&arrayFit);
353     delete his2D;
354     his1D = (TH1*) arrayFit.At(1);
355     his1D->Fit(fpol1,"QNRROB=0.8","QNR",-0.8,-0.1);
356     offsetdRC=fpol1->GetParameter(0);
357     slopedRC=fpol1->GetParameter(1);
358     offsetdRCchi2=fpol1->GetChisquare();
359     slopedRCchi2=fpol1->GetChisquare();
360     //
361     his1D->Fit(fpol1,"QNRROB=0.8","QNR",0.1,0.8);
362     offsetdRA=fpol1->GetParameter(0);
363     slopedRA=fpol1->GetParameter(1);
364     offsetdRAErr=fpol1->GetParError(0);
365     slopedRAErr=fpol1->GetParError(1);
366     offsetdRAchi2=fpol1->GetChisquare();
367     slopedRAchi2=fpol1->GetChisquare();
368     //
369     printf("DCA R QA report\n");
370     printf("offsetdRA\t%f\n",offsetdRA);
371     printf("slopedRA\t%f\n",slopedRA);
372     printf("offsetdRC\t%f\n",offsetdRC);
373     printf("slopedRC\t%f\n",slopedRC);
374     //
375     // dump values
376     //
377     (*pcstream)<<"tpcQA"<<
378         "offsetdRA="<< offsetdRA<<
379         "slopedRA="<< slopedRA<<
380         "offsetdRC="<< offsetdRC<<
381         "slopedRC="<<slopedRC<<
382         //
383         "offsetdRAErr="<< offsetdRAErr<<
384         "slopedRAErr="<< slopedRAErr<<
385         "offsetdRCErr="<< offsetdRCErr<<
386         "slopedRCErr="<<slopedRCErr<<
387         //
388         "offsetdRAchi2="<< offsetdRAchi2<<
389         "slopedRAchi2="<< slopedRAchi2<<
390         "offsetdRCchi2="<< offsetdRCchi2<<
391         "slopedRCchi2="<<slopedRCchi2;
392         
393     return 0;
394 }
395
396 //_____________________________________________________________________________
397 Int_t AliTPCPerformanceSummary::AnalyzeDCARPhiPos(const AliPerformanceTPC* pTPC, TTreeSRedirector* pcstream)
398 {
399     //
400     // Analyse DCA R imperfections for positive particles
401     //
402     
403     if (!pcstream) return 16;
404     if (!pTPC) return 16;
405
406     // variables:
407     static Double_t offsetdRAPos=0;
408     static Double_t slopedRAPos=0;
409     static Double_t offsetdRCPos=0;
410     static Double_t slopedRCPos=0;
411     static Double_t offsetdRAErrPos=0;
412     static Double_t slopedRAErrPos=0;
413     static Double_t offsetdRCErrPos=0;
414     static Double_t slopedRCErrPos=0;
415     static Double_t offsetdRAchi2Pos=0;
416     static Double_t slopedRAchi2Pos=0;
417     static Double_t offsetdRCchi2Pos=0;
418     static Double_t slopedRCchi2Pos=0;
419
420     //AliPerformanceTPC* pTPC =  dynamic_cast<AliPerformanceTPC*>(pTPCObject);    
421     TH1* his1D=0;
422     TH2* his2D=0;
423     static TF1 *fpol1 = new TF1("fpol1","pol1");
424     TObjArray arrayFit;
425     pTPC->GetTPCTrackHisto()->GetAxis(8)->SetRangeUser(0,10);
426     his2D = pTPC->GetTPCTrackHisto()->Projection(3,5); 
427     his2D->FitSlicesY(0,0,-1,10,"QNR",&arrayFit);
428     delete his2D;
429     pTPC->GetTPCTrackHisto()->GetAxis(8)->SetRangeUser(-10,10);
430     his1D = (TH1*) arrayFit.At(1);
431     his1D->Fit(fpol1,"QNRROB=0.8","QNR",-0.8,-0.1);
432     offsetdRCPos=fpol1->GetParameter(0);
433     slopedRCPos=fpol1->GetParameter(1);
434     offsetdRCchi2Pos=fpol1->GetChisquare();
435     slopedRCchi2Pos=fpol1->GetChisquare();
436     //
437     his1D->Fit(fpol1,"QNRROB=0.8","QNR",0.1,0.8);
438     offsetdRAPos=fpol1->GetParameter(0);
439     slopedRAPos=fpol1->GetParameter(1);
440     offsetdRAErrPos=fpol1->GetParError(0);
441     slopedRAErrPos=fpol1->GetParError(1);
442     offsetdRAchi2Pos=fpol1->GetChisquare();
443     slopedRAchi2Pos=fpol1->GetChisquare();
444     //
445     printf("DCA R QA Pos report\n");
446     printf("offsetdRAPos\t%f\n",offsetdRAPos);
447     printf("slopedRAPos\t%f\n",slopedRAPos);
448     printf("offsetdRCPos\t%f\n",offsetdRCPos);
449     printf("slopedRCPos\t%f\n",slopedRCPos);
450     //
451     // dump values
452     //
453     (*pcstream)<<"tpcQA"<<
454         "offsetdRAPos="<< offsetdRAPos<<
455         "slopedRAPos="<< slopedRAPos<<
456         "offsetdRCPos="<< offsetdRCPos<<
457         "slopedRCPos="<<slopedRCPos<<
458         //
459         "offsetdRAErrPos="<< offsetdRAErrPos<<
460         "slopedRAErrPos="<< slopedRAErrPos<<
461         "offsetdRCErrPos="<< offsetdRCErrPos<<
462         "slopedRCErrPos="<<slopedRCErrPos<<
463         //
464         "offsetdRAchi2Pos="<< offsetdRAchi2Pos<<
465         "slopedRAchi2Pos="<< slopedRAchi2Pos<<
466         "offsetdRCchi2Pos="<< offsetdRCchi2Pos<<
467         "slopedRCchi2Pos="<<slopedRCchi2Pos;
468         
469     return 0;
470 }
471
472 //_____________________________________________________________________________
473 Int_t AliTPCPerformanceSummary::AnalyzeDCARPhiNeg(const AliPerformanceTPC* pTPC, TTreeSRedirector* pcstream)
474 {
475     //
476     // Analyse DCA R imperfections for negative particles
477     //
478     if (!pcstream) return 32;
479     if (!pTPC) return 32;
480
481     // variables:
482     static Double_t offsetdRANeg=0;
483     static Double_t slopedRANeg=0;
484     static Double_t offsetdRCNeg=0;
485     static Double_t slopedRCNeg=0;
486     static Double_t offsetdRAErrNeg=0;
487     static Double_t slopedRAErrNeg=0;
488     static Double_t offsetdRCErrNeg=0;
489     static Double_t slopedRCErrNeg=0;
490     static Double_t offsetdRAchi2Neg=0;
491     static Double_t slopedRAchi2Neg=0;
492     static Double_t offsetdRCchi2Neg=0;
493     static Double_t slopedRCchi2Neg=0;
494
495     //AliPerformanceTPC* pTPC =  dynamic_cast<AliPerformanceTPC*>(pTPCObject);    
496     TH1* his1D=0;
497     TH2* his2D=0;
498     static TF1 *fpol1 = new TF1("fpol1","pol1");
499     TObjArray arrayFit;
500     pTPC->GetTPCTrackHisto()->GetAxis(8)->SetRangeUser(-10,0);
501     his2D = pTPC->GetTPCTrackHisto()->Projection(3,5); 
502     his2D->FitSlicesY(0,0,-1,10,"QNR",&arrayFit);
503     delete his2D;
504     pTPC->GetTPCTrackHisto()->GetAxis(8)->SetRangeUser(-10,10);
505     his1D = (TH1*) arrayFit.At(1);
506     his1D->Fit(fpol1,"QNRROB=0.8","QNR",-0.8,-0.1);
507     offsetdRCNeg=fpol1->GetParameter(0);
508     slopedRCNeg=fpol1->GetParameter(1);
509     offsetdRCchi2Neg=fpol1->GetChisquare();
510     slopedRCchi2Neg=fpol1->GetChisquare();
511     //
512     his1D->Fit(fpol1,"QNRROB=0.8","QNR",0.1,0.8);
513     offsetdRANeg=fpol1->GetParameter(0);
514     slopedRANeg=fpol1->GetParameter(1);
515     offsetdRAErrNeg=fpol1->GetParError(0);
516     slopedRAErrNeg=fpol1->GetParError(1);
517     offsetdRAchi2Neg=fpol1->GetChisquare();
518     slopedRAchi2Neg=fpol1->GetChisquare();
519     //
520     printf("DCA R QA Neg report\n");
521     printf("offsetdRANeg\t%f\n",offsetdRANeg);
522     printf("slopedRANeg\t%f\n",slopedRANeg);
523     printf("offsetdRCNeg\t%f\n",offsetdRCNeg);
524     printf("slopedRCNeg\t%f\n",slopedRCNeg);
525     //
526     // dump drift QA values
527     //
528     (*pcstream)<<"tpcQA"<<
529         "offsetdRANeg="<< offsetdRANeg<<
530         "slopedRANeg="<< slopedRANeg<<
531         "offsetdRCNeg="<< offsetdRCNeg<<
532         "slopedRCNeg="<<slopedRCNeg<<
533         //
534         "offsetdRAErrNeg="<< offsetdRAErrNeg<<
535         "slopedRAErrNeg="<< slopedRAErrNeg<<
536         "offsetdRCErrNeg="<< offsetdRCErrNeg<<
537         "slopedRCErrNeg="<<slopedRCErrNeg<<
538         //
539         "offsetdRAchi2Neg="<< offsetdRAchi2Neg<<
540         "slopedRAchi2Neg="<< slopedRAchi2Neg<<
541         "offsetdRCchi2Neg="<< offsetdRCchi2Neg<<
542         "slopedRCchi2Neg="<<slopedRCchi2Neg;
543         
544     return 0;
545 }
546
547 //_____________________________________________________________________________
548 Int_t AliTPCPerformanceSummary::AnalyzeNCL(const AliPerformanceTPC* pTPC, TTreeSRedirector* pcstream)
549 {
550     //
551     // Analyse number of TPC clusters 
552     //
553     
554     if (!pcstream) return 1;
555     if (!pTPC) return 1;
556  
557     // variables:
558     static Double_t meanTPCnclF=0;
559     static Double_t rmsTPCnclF=0;
560     static Double_t meanTPCChi2=0;
561     static Double_t rmsTPCChi2=0;  
562     static Double_t slopeATPCnclF=0;
563     static Double_t slopeCTPCnclF=0;
564     static Double_t slopeATPCnclFErr=0;
565     static Double_t slopeCTPCnclFErr=0;
566     static Double_t meanTPCncl=0;
567     static Double_t rmsTPCncl=0;
568     static Double_t slopeATPCncl=0;
569     static Double_t slopeCTPCncl=0;
570     static Double_t slopeATPCnclErr=0;
571     static Double_t slopeCTPCnclErr=0;  
572     TH1* his1D=0;
573     TProfile* hprof=0;
574     static TF1 *fpol1 = new TF1("fpol1","pol1");
575     //
576     // all clusters
577     // eta cut - +-1
578     // pt cut  - +-0.250 GeV
579     pTPC->GetTPCTrackHisto()->GetAxis(5)->SetRangeUser(-1.,1.);
580     pTPC->GetTPCTrackHisto()->GetAxis(7)->SetRangeUser(0.25,10);
581     his1D = pTPC->GetTPCTrackHisto()->Projection(0);
582     meanTPCncl= his1D->GetMean();
583     rmsTPCncl= his1D->GetRMS();
584     delete his1D;
585     his1D = pTPC->GetTPCTrackHisto()->Projection(1);
586     meanTPCChi2= his1D->GetMean();
587     rmsTPCChi2= his1D->GetRMS();
588     delete his1D;  
589     hprof = pTPC->GetTPCTrackHisto()->Projection(0,5)->ProfileX();
590     hprof->Fit(fpol1,"QNR","QNR",0,0.8);
591     slopeATPCncl= fpol1->GetParameter(1);
592     slopeATPCnclErr= fpol1->GetParError(1);
593     hprof->Fit(fpol1,"QNR","QNR",-0.8,0.0);
594     slopeCTPCncl= fpol1->GetParameter(1);
595     slopeCTPCnclErr= fpol1->GetParameter(1);
596     delete hprof;
597     //
598     // findable clusters
599     //
600     pTPC->GetTPCTrackHisto()->GetAxis(2)->SetRangeUser(0.4,1.1);
601     his1D = pTPC->GetTPCTrackHisto()->Projection(2);
602     meanTPCnclF= his1D->GetMean();
603     rmsTPCnclF= his1D->GetRMS();
604     delete his1D;
605     his1D = pTPC->GetTPCTrackHisto()->Projection(2,5)->ProfileX();
606     his1D->Fit(fpol1,"QNR","QNR",0,0.8);
607     slopeATPCnclF= fpol1->GetParameter(1);
608     slopeATPCnclFErr= fpol1->GetParError(1);
609     his1D->Fit(fpol1,"QNR","QNR",-0.8,0.0);
610     slopeCTPCnclF= fpol1->GetParameter(1);
611     slopeCTPCnclFErr= fpol1->GetParameter(1);
612     delete his1D;
613     printf("Cluster QA report\n");
614     printf("meanTPCnclF=\t%f\n",meanTPCnclF);
615     printf("rmsTPCnclF=\t%f\n",rmsTPCnclF);
616     printf("slopeATPCnclF=\t%f\n",slopeATPCnclF);
617     printf("slopeCTPCnclF=\t%f\n",slopeCTPCnclF);
618     printf("meanTPCncl=\t%f\n",meanTPCncl);
619     printf("rmsTPCncl=\t%f\n",rmsTPCncl);
620     printf("slopeATPCncl=\t%f\n",slopeATPCncl);
621     printf("slopeCTPCncl=\t%f\n",slopeCTPCncl);
622     printf("meanTPCChi2=\t%f\n",meanTPCChi2);
623     printf("rmsTPCChi2=\t%f\n",rmsTPCChi2);
624     //
625     // dump results to the tree
626     //
627     (*pcstream)<<"tpcQA"<<
628       "meanTPCnclF="<<meanTPCnclF <<   
629       "rmsTPCnclF="<<rmsTPCnclF <<
630       "meanTPCChi2="<<meanTPCChi2 <<
631       "rmsTPCChi2="<<rmsTPCChi2 <<
632       "slopeATPCnclF="<< slopeATPCnclF<<
633       "slopeCTPCnclF="<< slopeCTPCnclF<<
634       "slopeATPCnclFErr="<< slopeATPCnclFErr<<
635       "slopeCTPCnclFErr="<< slopeCTPCnclFErr<<
636       "meanTPCncl="<<meanTPCncl <<
637       "rmsTPCncl="<< rmsTPCncl<<
638       "slopeATPCncl="<< slopeATPCncl<<
639       "slopeCTPCncl="<< slopeCTPCncl<<
640       "slopeATPCnclErr="<< slopeATPCnclErr<<
641       "slopeCTPCnclErr="<< slopeCTPCnclErr;
642     
643     return 0;
644 }
645
646 //_____________________________________________________________________________
647 Int_t AliTPCPerformanceSummary::AnalyzeDrift(const AliPerformanceTPC* pTPC, TTreeSRedirector* pcstream)
648 {
649     //
650     // Analyse DCA Z imperferctions (drift velocity)
651     //
652     
653     if (!pcstream) return 2;
654     if (!pTPC) return 2;
655
656     // variables:
657     static Double_t offsetdZA=0;
658     static Double_t slopedZA=0;
659     static Double_t offsetdZC=0;
660     static Double_t slopedZC=0;
661     static Double_t offsetdZAErr=0;
662     static Double_t slopedZAErr=0;
663     static Double_t offsetdZCErr=0;
664     static Double_t slopedZCErr=0;
665     static Double_t offsetdZAchi2=0;
666     static Double_t slopedZAchi2=0;
667     static Double_t offsetdZCchi2=0;
668     static Double_t slopedZCchi2=0;
669     TH1* his1D=0;
670     TH2* his2D=0;
671     static TF1 *fpol1 = new TF1("fpol1","pol1");
672     TObjArray arrayFit;
673     his2D = pTPC->GetTPCTrackHisto()->Projection(4,5);
674     his2D->FitSlicesY(0,0,-1,10,"QNR",&arrayFit);
675     delete his2D;
676     his1D = (TH1*) arrayFit.At(1);
677     his1D->Fit(fpol1,"QNRROB=0.8","QNR",-0.8,-0.1);
678     offsetdZC=fpol1->GetParameter(0);
679     slopedZC=fpol1->GetParameter(1);
680     offsetdZCErr=fpol1->GetParError(0);
681     slopedZCErr=fpol1->GetParError(1);        
682     offsetdZCchi2=fpol1->GetChisquare();
683     slopedZCchi2=fpol1->GetChisquare();
684     //
685     his1D->Fit(fpol1,"QNRROB=0.8","QNR",0.1,0.8);
686     offsetdZA=fpol1->GetParameter(0);
687     slopedZA=fpol1->GetParameter(1);
688     offsetdZAErr=fpol1->GetParError(0);
689     slopedZAErr=fpol1->GetParError(1);
690     offsetdZAchi2=fpol1->GetChisquare();
691     slopedZAchi2=fpol1->GetChisquare();
692     //
693     printf("Drift velocity QA report\n");
694     printf("offsetdZA\t%f\n",offsetdZA);
695     printf("slopedZA\t%f\n",slopedZA);
696     printf("offsetdZC\t%f\n",offsetdZC);
697     printf("slopedZC\t%f\n",slopedZC);
698     //
699     // dump drift QA values
700     //
701     (*pcstream)<<"tpcQA"<<
702         "offsetdZA="<< offsetdZA<<
703         "slopedZA="<< slopedZA<<
704         "offsetdZC="<< offsetdZC<<
705         "slopedZC="<<slopedZC<<
706         //
707         "offsetdZAErr="<< offsetdZAErr<<
708         "slopedZAErr="<< slopedZAErr<<
709         "offsetdZCErr="<< offsetdZCErr<<
710         "slopedZCErr="<<slopedZCErr<<
711         //
712         "offsetdZAchi2="<< offsetdZAchi2<<
713         "slopedZAchi2="<< slopedZAchi2<<
714         "offsetdZCchi2="<< offsetdZCchi2<<
715         "slopedZCchi2="<<slopedZCchi2;
716     
717     return 0;
718 }
719
720 //_____________________________________________________________________________
721 Int_t AliTPCPerformanceSummary::AnalyzeDriftPos(const AliPerformanceTPC* pTPC, TTreeSRedirector* pcstream)
722 {
723     //
724     // Analyse DCA Z imperferctions (drift velocity)
725     // for positive particles
726     //
727     if (!pcstream) return 64;
728     if (!pTPC) return 64;
729         
730     // variables:
731     static Double_t offsetdZAPos=0;
732     static Double_t slopedZAPos=0;
733     static Double_t offsetdZCPos=0;
734     static Double_t slopedZCPos=0;
735     static Double_t offsetdZAErrPos=0;
736     static Double_t slopedZAErrPos=0;
737     static Double_t offsetdZCErrPos=0;
738     static Double_t slopedZCErrPos=0;
739     static Double_t offsetdZAchi2Pos=0;
740     static Double_t slopedZAchi2Pos=0;
741     static Double_t offsetdZCchi2Pos=0;
742     static Double_t slopedZCchi2Pos=0;
743     TH1* his1D=0;
744     TH2* his2D=0;
745     static TF1 *fpol1 = new TF1("fpol1","pol1");
746     TObjArray arrayFit;
747     pTPC->GetTPCTrackHisto()->GetAxis(8)->SetRangeUser(0,10);
748     his2D = pTPC->GetTPCTrackHisto()->Projection(4,5);
749     his2D->FitSlicesY(0,0,-1,10,"QNR",&arrayFit);
750     delete his2D;
751     pTPC->GetTPCTrackHisto()->GetAxis(8)->SetRangeUser(-10,10);
752     his1D = (TH1*) arrayFit.At(1);
753     his1D->Fit(fpol1,"QNRROB=0.8","QNR",-0.8,-0.1);
754     offsetdZCPos=fpol1->GetParameter(0);
755     slopedZCPos=fpol1->GetParameter(1);
756     offsetdZCErrPos=fpol1->GetParError(0);
757     slopedZCErrPos=fpol1->GetParError(1);        
758     offsetdZCchi2Pos=fpol1->GetChisquare();
759     slopedZCchi2Pos=fpol1->GetChisquare();
760     //
761     his1D->Fit(fpol1,"QNRROB=0.8","QNR",0.1,0.8);
762     offsetdZAPos=fpol1->GetParameter(0);
763     slopedZAPos=fpol1->GetParameter(1);
764     offsetdZAErrPos=fpol1->GetParError(0);
765     slopedZAErrPos=fpol1->GetParError(1);
766     offsetdZAchi2Pos=fpol1->GetChisquare();
767     slopedZAchi2Pos=fpol1->GetChisquare();
768     //
769     printf("Drift velocity QA report\n");
770     printf("offsetdZAPos\t%f\n",offsetdZAPos);
771     printf("slopedZAPos\t%f\n",slopedZAPos);
772     printf("offsetdZCPos\t%f\n",offsetdZCPos);
773     printf("slopedZCPos\t%f\n",slopedZCPos);
774     //
775     // dump drift QA values
776     //
777     (*pcstream)<<"tpcQA"<<
778         "offsetdZAPos="<< offsetdZAPos<<
779         "slopedZAPos="<< slopedZAPos<<
780         "offsetdZCPos="<< offsetdZCPos<<
781         "slopedZCPos="<<slopedZCPos<<
782         //
783         "offsetdZAErrPos="<< offsetdZAErrPos<<
784         "slopedZAErrPos="<< slopedZAErrPos<<
785         "offsetdZCErrPos="<< offsetdZCErrPos<<
786         "slopedZCErrPos="<<slopedZCErrPos<<
787         //
788         "offsetdZAchi2Pos="<< offsetdZAchi2Pos<<
789         "slopedZAchi2Pos="<< slopedZAchi2Pos<<
790         "offsetdZCchi2Pos="<< offsetdZCchi2Pos<<
791         "slopedZCchi2Pos="<<slopedZCchi2Pos;
792         
793     return 0;
794 }
795
796 //_____________________________________________________________________________
797 Int_t AliTPCPerformanceSummary::AnalyzeDriftNeg(const AliPerformanceTPC* pTPC, TTreeSRedirector* pcstream)
798 {
799     //
800     // Analyse DCA Z imperferctions (drift velocity)
801     // for negative particles
802     //
803     if (!pcstream) return 128;
804     if (!pTPC) return 128;
805             
806     // variables:
807     static Double_t offsetdZANeg=0;
808     static Double_t slopedZANeg=0;
809     static Double_t offsetdZCNeg=0;
810     static Double_t slopedZCNeg=0;
811     static Double_t offsetdZAErrNeg=0;
812     static Double_t slopedZAErrNeg=0;
813     static Double_t offsetdZCErrNeg=0;
814     static Double_t slopedZCErrNeg=0;
815     static Double_t offsetdZAchi2Neg=0;
816     static Double_t slopedZAchi2Neg=0;
817     static Double_t offsetdZCchi2Neg=0;
818     static Double_t slopedZCchi2Neg=0;
819     TH1* his1D=0;
820     TH2* his2D=0;
821     static TF1 *fpol1 = new TF1("fpol1","pol1");
822     TObjArray arrayFit;
823     pTPC->GetTPCTrackHisto()->GetAxis(8)->SetRangeUser(-10,0);
824     his2D = pTPC->GetTPCTrackHisto()->Projection(4,5);
825     his2D->FitSlicesY(0,0,-1,10,"QNR",&arrayFit);
826     delete his2D;
827     pTPC->GetTPCTrackHisto()->GetAxis(8)->SetRangeUser(-10,10);
828     his1D = (TH1*) arrayFit.At(1);
829     his1D->Fit(fpol1,"QNRROB=0.8","QNR",-0.8,-0.1);
830     offsetdZCNeg=fpol1->GetParameter(0);
831     slopedZCNeg=fpol1->GetParameter(1);
832     offsetdZCErrNeg=fpol1->GetParError(0);
833     slopedZCErrNeg=fpol1->GetParError(1);        
834     offsetdZCchi2Neg=fpol1->GetChisquare();
835     slopedZCchi2Neg=fpol1->GetChisquare();
836     //
837     his1D->Fit(fpol1,"QNRROB=0.8","QNR",0.1,0.8);
838     offsetdZANeg=fpol1->GetParameter(0);
839     slopedZANeg=fpol1->GetParameter(1);
840     offsetdZAErrNeg=fpol1->GetParError(0);
841     slopedZAErrNeg=fpol1->GetParError(1);
842     offsetdZAchi2Neg=fpol1->GetChisquare();
843     slopedZAchi2Neg=fpol1->GetChisquare();
844     //
845     printf("Drift velocity QA report\n");
846     printf("offsetdZANeg\t%f\n",offsetdZANeg);
847     printf("slopedZANeg\t%f\n",slopedZANeg);
848     printf("offsetdZCNeg\t%f\n",offsetdZCNeg);
849     printf("slopedZCNeg\t%f\n",slopedZCNeg);
850     //
851     // dump drift QA values
852     //
853     (*pcstream)<<"tpcQA"<<
854         "offsetdZANeg="<< offsetdZANeg<<
855         "slopedZANeg="<< slopedZANeg<<
856         "offsetdZCNeg="<< offsetdZCNeg<<
857         "slopedZCNeg="<<slopedZCNeg<<
858         //
859         "offsetdZAErrNeg="<< offsetdZAErrNeg<<
860         "slopedZAErrNeg="<< slopedZAErrNeg<<
861         "offsetdZCErrNeg="<< offsetdZCErrNeg<<
862         "slopedZCErrNeg="<<slopedZCErrNeg<<
863         //
864         "offsetdZAchi2Neg="<< offsetdZAchi2Neg<<
865         "slopedZAchi2Neg="<< slopedZAchi2Neg<<
866         "offsetdZCchi2Neg="<< offsetdZCchi2Neg<<
867         "slopedZCchi2Neg="<<slopedZCchi2Neg;
868     
869     return 0;
870 }
871
872 //_____________________________________________________________________________
873 Int_t AliTPCPerformanceSummary::AnalyzeGain(const AliPerformanceDEdx* pTPCgain, TTreeSRedirector* pcstream)
874 {
875     //
876     // Analyse Gain
877     //
878     
879     if (!pcstream) return 4;
880     if (!pTPCgain) return 4;
881
882     static TVectorD meanMIPvsSector(36);
883     static TVectorD sector(36);
884     static Float_t meanMIP = 0;
885     static Float_t resolutionMIP = 0;
886     static Float_t attachSlopeC = 0;
887     static Float_t attachSlopeA = 0;
888
889     meanMIPvsSector.Zero();
890     //
891     // select MIP particles
892     //
893     pTPCgain->GetDeDxHisto()->GetAxis(7)->SetRangeUser(0.4,0.55);
894     pTPCgain->GetDeDxHisto()->GetAxis(0)->SetRangeUser(35,60);
895     pTPCgain->GetDeDxHisto()->GetAxis(6)->SetRangeUser(80,160);
896     pTPCgain->GetDeDxHisto()->GetAxis(5)->SetRangeUser(-1,1);
897     //
898     // MIP position and resolution
899     //
900     TF1 gausFit("gausFit","gaus");
901     TH1 * hisProj1D = pTPCgain->GetDeDxHisto()->Projection(0);
902     hisProj1D->Fit(&gausFit,"QN","QN");
903     delete hisProj1D;
904     meanMIP = gausFit.GetParameter(1);
905     resolutionMIP = 0;
906     if (meanMIP!=0) resolutionMIP = gausFit.GetParameter(2)/meanMIP;
907     //
908     // MIP position vs. dip angle (attachment)
909     //
910     pTPCgain->GetDeDxHisto()->GetAxis(5)->SetRangeUser(-3,0);
911     TH2* his2D = pTPCgain->GetDeDxHisto()->Projection(0,5);
912     TF1 * fpol = new TF1("fpol","pol1");
913     TObjArray arrayFit;
914     his2D->FitSlicesY(0,0,-1,10,"QN",&arrayFit);
915     delete his2D;
916     TH1 * his1D = (TH1*) arrayFit.At(1);
917     his1D->Fit(fpol,"QNROB=0.8","QNR",-1,0);
918     attachSlopeC = fpol->GetParameter(1);
919     //
920     pTPCgain->GetDeDxHisto()->GetAxis(5)->SetRangeUser(0,3);
921     TH2* his2DA = pTPCgain->GetDeDxHisto()->Projection(0,5);
922     TF1 * fpolA = new TF1("fpolA","pol1");
923     TObjArray arrayFitA;
924     his2DA->FitSlicesY(0,0,-1,10,"QN",&arrayFit);
925     delete his2DA;
926     TH1 * his1DA = (TH1*) arrayFit.At(1);
927     his1DA->Fit(fpolA,"QNROB=0.8","QN",0,1);
928     attachSlopeA = fpolA->GetParameter(1);
929     //
930     // MIP position vs. sector
931     //
932     pTPCgain->GetDeDxHisto()->GetAxis(5)->SetRangeUser(-3,0); // C side
933     for(Int_t i = 0; i < 18; i++) { // loop over sectors; correct mapping to be checked!
934         //TH1* his1D=0;
935         Float_t phiLow = -TMath::Pi() + i*(20./360.)*(2*TMath::Pi());
936         Float_t phiUp    = -TMath::Pi() + (i+1)*(20./360.)*(2*TMath::Pi());
937         pTPCgain->GetDeDxHisto()->GetAxis(1)->SetRangeUser(phiLow,phiUp);
938         his1D = pTPCgain->GetDeDxHisto()->Projection(0);
939         TF1 gausFunc("gausFunc","gaus");
940         his1D->Fit(&gausFunc, "QN");
941         meanMIPvsSector(i) = gausFunc.GetParameter(1);
942         sector(i)=i;
943         delete his1D;
944     }
945     //
946     pTPCgain->GetDeDxHisto()->GetAxis(5)->SetRangeUser(0,3); // A side
947     for(Int_t i = 0; i < 18; i++) { // loop over sectors; correct mapping to be checked!
948         TH1* his1D=0;
949         Float_t phiLow = -TMath::Pi() + i*(20./360.)*(2*TMath::Pi());
950         Float_t phiUp    = -TMath::Pi() + (i+1)*(20./360.)*(2*TMath::Pi());
951         pTPCgain->GetDeDxHisto()->GetAxis(1)->SetRangeUser(phiLow,phiUp);
952         his1D = pTPCgain->GetDeDxHisto()->Projection(0);
953         TF1 gausFunc("gausFunc","gaus");
954         his1D->Fit(&gausFunc, "QN");
955         meanMIPvsSector(i+18) = gausFunc.GetParameter(1);
956         sector(i+18)=i+18;
957         delete his1D;
958     }
959     //
960     printf("Gain QA report\n");
961     printf("MIP mean\t%f\n",meanMIP);
962     printf("MIP resolution\t%f\n",resolutionMIP);
963     printf("MIPslopeA\t%f\n",attachSlopeA);
964     printf("MIPslopeC\t%f\n",attachSlopeC);
965     //
966     (*pcstream)<<"tpcQA"<<
967         "MIPattachSlopeC="<<attachSlopeC<<
968         "MIPattachSlopeA="<<attachSlopeA<<
969         "resolutionMIP="<<resolutionMIP<<
970         "meanMIPvsSector.="<<&meanMIPvsSector<<
971         "sector.="<<&sector<<
972         "meanMIP="<<meanMIP;
973
974     return 0;
975 }
976
977 //_____________________________________________________________________________
978 Int_t AliTPCPerformanceSummary::AnalyzeEvent(const AliPerformanceTPC* pTPC, TTreeSRedirector* pcstream)
979 {
980     //
981     // Analyse Primary Vertex Distribution and Multiplicities
982     //
983   if (!pcstream) return 1;
984     if (!pTPC) return 1;
985     //
986     // 
987     //
988     static Double_t meanVertX=0;
989     static Double_t rmsVertX=0;
990     static Double_t meanVertY=0;
991     static Double_t rmsVertY=0;
992     static Double_t meanVertZ=0;
993     static Double_t rmsVertZ=0;
994     static Double_t vertStatus=0;
995     static Double_t meanMult=0;
996     static Double_t rmsMult=0;
997     static Double_t meanMultPos=0;
998     static Double_t rmsMultPos=0;
999     static Double_t meanMultNeg=0;
1000     static Double_t rmsMultNeg=0;
1001     static Double_t vertAll = 0;
1002     static Double_t vertOK = 0;
1003     
1004     TH1* his1D=0;
1005     
1006     his1D = pTPC->GetTPCEventHisto()->Projection(6);
1007     vertAll = his1D->GetEntries();
1008     delete his1D;
1009     
1010     pTPC->GetTPCEventHisto()->GetAxis(6)->SetRangeUser(1,1);
1011     
1012     his1D = pTPC->GetTPCEventHisto()->Projection(6);
1013     vertOK = his1D->GetEntries();
1014     delete his1D;
1015     if (vertAll>=1) {
1016             vertStatus = vertOK / vertAll;
1017     }
1018     his1D = pTPC->GetTPCEventHisto()->Projection(0);
1019     meanVertX = his1D->GetMean();    
1020     rmsVertX    = his1D->GetRMS();
1021     delete his1D;
1022     
1023     his1D = pTPC->GetTPCEventHisto()->Projection(1);
1024     meanVertY = his1D->GetMean();
1025     rmsVertY    = his1D->GetRMS();
1026     delete his1D;
1027     
1028     his1D = pTPC->GetTPCEventHisto()->Projection(2);
1029     meanVertZ = his1D->GetMean();
1030     rmsVertZ    = his1D->GetRMS();
1031     delete his1D;
1032     
1033     his1D = pTPC->GetTPCEventHisto()->Projection(3);
1034     meanMult    = his1D->GetMean();
1035     rmsMult     = his1D->GetRMS();
1036     delete his1D;
1037     
1038     his1D = pTPC->GetTPCEventHisto()->Projection(4);
1039     meanMultPos    = his1D->GetMean();
1040     rmsMultPos     = his1D->GetRMS();
1041     delete his1D;
1042     
1043     his1D = pTPC->GetTPCEventHisto()->Projection(5);
1044     meanMultNeg    = his1D->GetMean();
1045     rmsMultNeg     = his1D->GetRMS();
1046     delete his1D;
1047     
1048     pTPC->GetTPCEventHisto()->GetAxis(6)->SetRangeUser(0,1);
1049     //
1050     (*pcstream)<<"tpcQA"<<
1051         "meanVertX="<<meanVertX<<
1052         "rmsVertX="<<rmsVertX<<
1053         "meanVertY="<<meanVertY<<
1054         "rmsVertY="<<rmsVertY<<
1055         "meanVertZ="<<meanVertZ<<
1056         "rmsVertZ="<<rmsVertZ<<
1057         "vertStatus="<<vertStatus<<
1058         "vertAll="<<vertAll<<
1059         "vertOK="<<vertOK<<
1060         "meanMult="<<meanMult<<
1061         "rmsMult="<<rmsMult<<
1062         "meanMultPos="<<meanMultPos<<
1063         "rmsMultPos="<<rmsMultPos<<
1064         "meanMultNeg="<<meanMultNeg<<
1065         "rmsMultNeg="<<rmsMultNeg;     
1066      
1067     return 0;
1068 }