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.
11 // by M.Knichel 15/10/2010
12 //------------------------------------------------------------------------------
27 #include "THnSparse.h"
34 #include "AliGRPObject.h"
35 #include "AliTPCcalibDB.h"
36 #include "AliTPCcalibDButil.h"
37 #include "TTreeStream.h"
38 #include "AliPerformanceTPC.h"
39 #include "AliPerformanceDEdx.h"
40 #include "AliPerformanceDCA.h"
41 #include "AliPerformanceMatch.h"
43 #include "AliTPCPerformanceSummary.h"
47 ClassImp(AliTPCPerformanceSummary)
49 Bool_t AliTPCPerformanceSummary::fgForceTHnSparse = kFALSE;
52 //_____________________________________________________________________________
53 void AliTPCPerformanceSummary::WriteToTTreeSRedirector(const AliPerformanceTPC* pTPC, const AliPerformanceDEdx* pTPCgain, const AliPerformanceMatch* pTPCMatch,const AliPerformanceMatch* pTPCPull, const AliPerformanceMatch* pConstrain, TTreeSRedirector* const pcstream, Int_t run)
56 // Extracts performance parameters from pTPC and pTPCgain.
57 // Output is written to pcstream.
58 // The run number must be provided since it is not stored in
59 // AliPerformanceTPC or AliPerformanceDEdx.
62 if (pTPCMatch) {run = pTPCMatch->GetRunNumber(); }
63 if (pTPCgain) {run = pTPCgain->GetRunNumber(); }
64 if (pTPC) { run = pTPC->GetRunNumber(); }
68 //AliTPCcalibDB *calibDB=0;
70 // AliTPCcalibDButil *dbutil =0;
76 //Float_t currentL3 =0;
77 //Int_t polarityL3 = 0;
80 //calibDB = AliTPCcalibDB::Instance();
82 // dbutil= new AliTPCcalibDButil;
84 //printf("Processing run %d ...\n",run);
86 //AliTPCcalibDB::Instance()->SetRun(run);
88 // dbutil->UpdateFromCalibDB();
89 // dbutil->SetReferenceRun(run);
90 // dbutil->UpdateRefDataFromOCDB();
92 //if (calibDB->GetGRP(run)){
93 //startTimeGRP = AliTPCcalibDB::GetGRP(run)->GetTimeStart();
94 //stopTimeGRP = AliTPCcalibDB::GetGRP(run)->GetTimeEnd();
95 //currentL3 = AliTPCcalibDB::GetL3Current(run);
96 //polarityL3 = AliTPCcalibDB::GetL3Polarity(run);
97 //bz = AliTPCcalibDB::GetBz(run);
100 //runType = AliTPCcalibDB::GetRunType(run).Data();
102 time = (startTimeGRP+stopTimeGRP)/2;
103 duration = (stopTimeGRP-startTimeGRP);
105 if (!pcstream) return;
106 (*pcstream)<<"tpcQA"<<
109 "startTimeGRP="<<startTimeGRP<<
110 "stopTimeGRP="<<stopTimeGRP<<
111 "duration="<<duration<<
112 "runType.="<<&runType;
114 pTPC->GetTPCTrackHisto()->GetAxis(9)->SetRangeUser(0.5,1.5);
115 pTPC->GetTPCTrackHisto()->GetAxis(7)->SetRangeUser(0.25,10);
116 pTPC->GetTPCTrackHisto()->GetAxis(5)->SetRangeUser(-1,1);
117 AnalyzeNCL(pTPC, pcstream);
118 AnalyzeDrift(pTPC, pcstream);
119 AnalyzeDriftPos(pTPC, pcstream);
120 AnalyzeDriftNeg(pTPC, pcstream);
121 AnalyzeDCARPhi(pTPC, pcstream);
122 AnalyzeDCARPhiPos(pTPC, pcstream);
123 AnalyzeDCARPhiNeg(pTPC, pcstream);
124 AnalyzeEvent(pTPC, pcstream);
126 AnalyzePt(pTPC,pcstream);
127 AnalyzeChargeOverPt(pTPC,pcstream);
128 AnalyzeQAPosNegDpT(pTPC,pcstream);
129 AnalyzeQADCAFitParameter(pTPC,pcstream);
131 pTPC->GetTPCTrackHisto()->GetAxis(9)->SetRangeUser(-10,10);
132 pTPC->GetTPCTrackHisto()->GetAxis(7)->SetRangeUser(0,100);
133 pTPC->GetTPCTrackHisto()->GetAxis(5)->SetRangeUser(-10,10);
135 AnalyzeGain(pTPCgain, pcstream);
136 AnalyzeMatch(pTPCMatch, pcstream);
137 AnalyzePull(pTPCPull, pcstream);
138 AnalyzeConstrain(pConstrain, pcstream);
140 (*pcstream)<<"tpcQA"<<"\n";
143 //_____________________________________________________________________________
144 void AliTPCPerformanceSummary::WriteToFile(const AliPerformanceTPC* pTPC, const AliPerformanceDEdx* pTPCgain, const AliPerformanceMatch* pMatch, const AliPerformanceMatch* pPull, const AliPerformanceMatch* pConstrain, const Char_t* outfile, Int_t run)
147 // Extracts performance parameters from pTPC and pTPCgain.
148 // Output is written to a TTree saved in outfile.
149 // The run number must be provided since it is not stored in
150 // AliPerformanceTPC or AliPerformanceDEdx.
152 // The function creates a TTreeSRedirector and calls the
153 // function WriteToTTreeSRedirector.
156 if (!outfile) return;
157 TTreeSRedirector* pcstream = 0;
158 pcstream = new TTreeSRedirector(outfile);
159 if (!pcstream) return;
160 WriteToTTreeSRedirector(pTPC, pTPCgain, pMatch, pPull, pConstrain, pcstream, run);
161 if (pcstream) { delete pcstream; pcstream = 0; }
165 //_____________________________________________________________________________
166 Int_t AliTPCPerformanceSummary::MakeReport(const Char_t* infile, const Char_t* outfile, Int_t run)
169 // Reads QA information (AliPerformanceTPC and AliPerformanceDEdx) from
170 // infile (this must be a rootfile) and writes the output to a TTree
171 // stored in outfile.
172 // The run number must be provided since it is not stored in
173 // AliPerformanceTPC or AliPerformanceDEdx.
175 // The input objects must be named "AliPerformanceTPC" and
176 // "AliPerformanceDEdxTPCInner" and stored in a TList which name must
177 // be one of the following: "TPC", "TPCQA", "TPC_PerformanceQA"
178 // or "TPC_PerformanceQA/TPC" (with directory)
181 if (!infile) return -1;
182 if (!outfile) return -1;
184 f=TFile::Open(infile,"read");
186 printf("File %s not available\n", infile);
190 list = dynamic_cast<TList*>(f->Get("TPC"));
191 if (!list) { list = dynamic_cast<TList*>(f->Get("TPCQA")); }
192 if (!list) { list = dynamic_cast<TList*>(f->Get("TPC_PerformanceQA/TPCQA")); }
193 if (!list) { list = dynamic_cast<TList*>(f->Get("TPC_PerformanceQA")); }
194 if (!list) { list = dynamic_cast<TList*>(f->Get("ITSTPCMatch")); }
196 printf("QA %s not available\n", infile);
199 AliPerformanceTPC* pTPC = 0;
200 AliPerformanceDEdx* pTPCgain = 0;
201 AliPerformanceMatch* pTPCmatch = 0;
202 AliPerformanceMatch* pTPCPull = 0;
203 AliPerformanceMatch* pConstrain = 0;
205 if (list) { pTPC = dynamic_cast<AliPerformanceTPC*>(list->FindObject("AliPerformanceTPC")); }
206 if (list) { pTPCgain = dynamic_cast<AliPerformanceDEdx*>(list->FindObject("AliPerformanceDEdxTPCInner")); }
207 if (list) { pTPCmatch = dynamic_cast<AliPerformanceMatch*>(list->FindObject("AliPerformanceMatchTPCITS")); }
208 if (list) { pTPCPull = dynamic_cast<AliPerformanceMatch*>(list->FindObject("AliPerformanceMatchITSTPC")); }
209 if (list) { pConstrain = dynamic_cast<AliPerformanceMatch*>(list->FindObject("AliPerformanceMatchTPCConstrain")); }
211 Int_t returncode = 0;
212 WriteToFile(pTPC, pTPCgain, pTPCmatch , pTPCPull, pConstrain, outfile, run);
213 if (f) { delete f; f=0; }
217 //_____________________________________________________________________________
218 Int_t AliTPCPerformanceSummary::ProduceTrends(const Char_t* infilelist, const Char_t* outfile)
221 // Produces trend graphs.
223 // Input: infilelist is a textfile with one rootfile per line.
224 // There should be one rootfile for each run, the rootfile must
225 // contain the output of the MakeReport function
226 // Output: the information for all runs is merged into a TTree
227 // that is saved in outfile along with the trend graphs.
228 // Trend graphs are stored as TCanvas objects to include axis labels etc.
231 if (!infilelist) return -1;
232 if (!outfile) return -1;
234 TChain* chain = new TChain("tpcQA");
235 if(!chain) return -1;
244 if (!currentFile.Contains("root")) continue; // protection
245 chain->Add(currentFile.Data());
248 //TTree *tree = chain;
249 TTree *tree = chain->CopyTree("1");
251 if (chain) { delete chain; chain=0; }
252 //TGraph* graph = dynamic_cast<TGraph*>(tree->DrawClone("run:run"));
253 //TGraph *graph = (TGraph*)gPad->GetPrimitive("Graph");
255 TFile* out = new TFile(outfile,"RECREATE");
259 const Char_t* condition = "meanTPCncl>0";
260 SaveGraph(tree,"meanTPCnclF","run",condition);
261 SaveGraph(tree,"rmsTPCnclF","run",condition);
262 SaveGraph(tree,"meanTPCChi2","run",condition);
263 SaveGraph(tree,"rmsTPCChi2","run",condition);
264 SaveGraph(tree,"slopeATPCnclF","run",condition);
265 SaveGraph(tree,"slopeCTPCnclF","run",condition);
266 SaveGraph(tree,"slopeATPCnclFErr","run",condition);
267 SaveGraph(tree,"slopeCTPCnclFErr","run",condition);
268 SaveGraph(tree,"meanTPCncl","run",condition);
269 SaveGraph(tree,"rmsTPCncl","run",condition);
270 SaveGraph(tree,"slopeATPCncl","run",condition);
271 SaveGraph(tree,"slopeCTPCncl","run",condition);
272 SaveGraph(tree,"slopeATPCnclErr","run",condition);
273 SaveGraph(tree,"slopeCTPCnclErr","run",condition);
275 SaveGraph(tree,"offsetdRA","run",condition);
276 SaveGraph(tree,"slopedRA","run",condition);
277 SaveGraph(tree,"offsetdRC","run",condition);
278 SaveGraph(tree,"slopedRC","run",condition);
279 SaveGraph(tree,"offsetdRAErr","run",condition);
280 SaveGraph(tree,"slopedRAErr","run",condition);
281 SaveGraph(tree,"offsetdRCErr","run",condition);
282 SaveGraph(tree,"slopedRCErr","run",condition);
283 SaveGraph(tree,"offsetdRAchi2","run",condition);
284 SaveGraph(tree,"slopedRAchi2","run",condition);
285 SaveGraph(tree,"offsetdRCchi2","run",condition);
286 SaveGraph(tree,"slopedRCchi2","run",condition);
288 SaveGraph(tree,"offsetdRAPos","run",condition);
290 SaveGraph(tree,"slopedRAPos","run",condition);
291 SaveGraph(tree,"offsetdRCPos","run",condition);
292 SaveGraph(tree,"slopedRCPos","run",condition);
293 SaveGraph(tree,"offsetdRAErrPos","run",condition);
294 SaveGraph(tree,"slopedRAErrPos","run",condition);
295 SaveGraph(tree,"offsetdRCErrPos","run",condition);
296 SaveGraph(tree,"slopedRCErrPos","run",condition);
297 SaveGraph(tree,"offsetdRAchi2Pos","run",condition);
298 SaveGraph(tree,"slopedRAchi2Pos","run",condition);
299 SaveGraph(tree,"offsetdRCchi2Pos","run",condition);
300 SaveGraph(tree,"slopedRCchi2Pos","run",condition);
302 SaveGraph(tree,"offsetdRANeg","run",condition);
303 SaveGraph(tree,"slopedRANeg","run",condition);
304 SaveGraph(tree,"offsetdRCNeg","run",condition);
305 SaveGraph(tree,"slopedRCNeg","run",condition);
306 SaveGraph(tree,"offsetdRAErrNeg","run",condition);
307 SaveGraph(tree,"slopedRAErrNeg","run",condition);
308 SaveGraph(tree,"offsetdRCErrNeg","run",condition);
309 SaveGraph(tree,"slopedRCErrNeg","run",condition);
310 SaveGraph(tree,"offsetdRAchi2Neg","run",condition);
311 SaveGraph(tree,"slopedRAchi2Neg","run",condition);
312 SaveGraph(tree,"offsetdRCchi2Neg","run",condition);
313 SaveGraph(tree,"slopedRCchi2Neg","run",condition);
315 SaveGraph(tree,"offsetdZAPos","run",condition);
316 SaveGraph(tree,"slopedZAPos","run",condition);
317 SaveGraph(tree,"offsetdZCPos","run",condition);
318 SaveGraph(tree,"slopedZCPos","run",condition);
319 SaveGraph(tree,"offsetdZAErrPos","run",condition);
320 SaveGraph(tree,"slopedZAErrPos","run",condition);
321 SaveGraph(tree,"offsetdZCErrPos","run",condition);
322 SaveGraph(tree,"slopedZCErrPos","run",condition);
323 SaveGraph(tree,"offsetdZAchi2Pos","run",condition);
324 SaveGraph(tree,"slopedZAchi2Pos","run",condition);
325 SaveGraph(tree,"offsetdZCchi2Pos","run",condition);
326 SaveGraph(tree,"slopedZCchi2Pos","run",condition);
328 SaveGraph(tree,"offsetdZANeg","run",condition);
329 SaveGraph(tree,"slopedZANeg","run",condition);
330 SaveGraph(tree,"offsetdZCNeg","run",condition);
331 SaveGraph(tree,"slopedZCNeg","run",condition);
332 SaveGraph(tree,"offsetdZAErrNeg","run",condition);
333 SaveGraph(tree,"slopedZAErrNeg","run",condition);
334 SaveGraph(tree,"offsetdZCErrNeg","run",condition);
335 SaveGraph(tree,"slopedZCErrNeg","run",condition);
336 SaveGraph(tree,"offsetdZAchi2Neg","run",condition);
337 SaveGraph(tree,"slopedZAchi2Neg","run",condition);
338 SaveGraph(tree,"offsetdZCchi2Neg","run",condition);
339 SaveGraph(tree,"slopedZCchi2Neg","run",condition);
341 SaveGraph(tree,"offsetdZA","run",condition);
342 SaveGraph(tree,"slopedZA","run",condition);
343 SaveGraph(tree,"offsetdZC","run",condition);
344 SaveGraph(tree,"slopedZC","run",condition);
345 SaveGraph(tree,"offsetdZAErr","run",condition);
346 SaveGraph(tree,"slopedZAErr","run",condition);
347 SaveGraph(tree,"offsetdZCErr","run",condition);
348 SaveGraph(tree,"slopedZCErr","run",condition);
349 SaveGraph(tree,"offsetdZAchi2","run",condition);
350 SaveGraph(tree,"slopedZAchi2","run",condition);
351 SaveGraph(tree,"offsetdZCchi2","run",condition);
352 SaveGraph(tree,"slopedZCchi2","run",condition);
354 SaveGraph(tree,"meanVertX","run",condition);
355 SaveGraph(tree,"rmsVertX","run",condition);
356 SaveGraph(tree,"meanVertY","run",condition);
357 SaveGraph(tree,"rmsVertY","run",condition);
358 SaveGraph(tree,"meanVertZ","run",condition);
359 SaveGraph(tree,"rmsVertZ","run",condition);
360 SaveGraph(tree,"vertStatus","run",condition);
361 SaveGraph(tree,"meanMult","run",condition);
362 SaveGraph(tree,"rmsMult","run",condition);
363 SaveGraph(tree,"meanMultPos","run",condition);
364 SaveGraph(tree,"rmsMultPos","run",condition);
365 SaveGraph(tree,"meanMultNeg","run",condition);
366 SaveGraph(tree,"rmsMultNeg","run",condition);
367 SaveGraph(tree,"vertAll","run",condition);
368 SaveGraph(tree,"vertOK","run",condition);
371 SaveGraph(tree,"meanPtAPos","run",condition);
372 SaveGraph(tree,"mediumPtAPos","run",condition);
373 SaveGraph(tree,"highPtAPos","run",condition);
374 SaveGraph(tree,"meanPtCPos","run",condition);
375 SaveGraph(tree,"mediumPtCPos","run",condition);
376 SaveGraph(tree,"highPtCPos","run",condition);
377 SaveGraph(tree,"meanPtANeg","run",condition);
378 SaveGraph(tree,"mediumPtANeg","run",condition);
379 SaveGraph(tree,"highPtANeg","run",condition);
380 SaveGraph(tree,"meanPtCNeg","run",condition);
381 SaveGraph(tree,"mediumPtCNeg","run",condition);
382 SaveGraph(tree,"highPtCNeg","run",condition);
384 SaveGraph(tree,"qOverPt","run",condition);
385 SaveGraph(tree,"qOverPtA","run",condition);
386 SaveGraph(tree,"qOverPtC","run",condition);
388 SaveGraph(tree,"dcarAP0","run",condition);
389 SaveGraph(tree,"dcarAP1","run",condition);
390 SaveGraph(tree,"dcarCP0","run",condition);
391 SaveGraph(tree,"dcarCP1","run",condition);
394 SaveGraph(tree,"tpcItsMatchA","run",condition);
395 SaveGraph(tree,"tpcItsMatchHighPtA","run",condition);
396 SaveGraph(tree,"tpcItsMatchC","run",condition);
397 SaveGraph(tree,"tpcItsMatchHighPtC","run",condition);
399 SaveGraph(tree,"phiPull","run",condition);
400 SaveGraph(tree,"phiPullHighPt","run",condition);
401 SaveGraph(tree,"ptPull","run",condition);
402 SaveGraph(tree,"ptPullHighPt","run",condition);
403 SaveGraph(tree,"yPull","run",condition);
404 SaveGraph(tree,"yPullHighPt","run",condition);
405 SaveGraph(tree,"zPull","run",condition);
406 SaveGraph(tree,"zPullHighPt","run",condition);
407 SaveGraph(tree,"lambdaPull","run",condition);
408 SaveGraph(tree,"lambdaPullHighPt","run",condition);
410 SaveGraph(tree,"tpcConstrainPhiA","run",condition);
411 SaveGraph(tree,"tpcConstrainPhiC","run",condition);
413 SaveGraph(tree,"deltaPt","run",condition);
414 SaveGraph(tree,"deltaPtchi2","run",condition);
415 SaveGraph(tree,"deltaPtA","run",condition);
416 SaveGraph(tree,"deltaPtchi2A","run",condition);
417 SaveGraph(tree,"deltaPtC","run",condition);
418 SaveGraph(tree,"deltaPtchi2C","run",condition);
419 SaveGraph(tree,"deltaPtA_Err","run",condition);
420 SaveGraph(tree,"deltaPtA_Err","run",condition);
421 SaveGraph(tree,"deltaPtC_Err","run",condition);
423 ////////////////////////////////////////////////////////////////////////////////////////////////////////
424 //save dca fit parameters
425 SaveGraph(tree,"dcar_posA_0","run",condition);
426 SaveGraph(tree,"dcar_posA_1","run",condition);
427 SaveGraph(tree,"dcar_posA_2","run",condition);
428 SaveGraph(tree,"dcar_posA_chi2","run",condition);
429 SaveGraph(tree,"dcar_posA_0_Err","run",condition);
430 SaveGraph(tree,"dcar_posA_1_Err","run",condition);
431 SaveGraph(tree,"dcar_posA_2_Err","run",condition);
433 SaveGraph(tree,"dcaz_posA_0","run",condition);
434 SaveGraph(tree,"dcaz_posA_1","run",condition);
435 SaveGraph(tree,"dcaz_posA_2","run",condition);
436 SaveGraph(tree,"dcaz_posA_chi2","run",condition);
437 SaveGraph(tree,"dcaz_posA_0_Err","run",condition);
438 SaveGraph(tree,"dcaz_posA_1_Err","run",condition);
439 SaveGraph(tree,"dcaz_posA_2_Err","run",condition);
441 SaveGraph(tree,"dcaz_posC_0","run",condition);
442 SaveGraph(tree,"dcaz_posC_1","run",condition);
443 SaveGraph(tree,"dcaz_posC_2","run",condition);
444 SaveGraph(tree,"dcaz_posC_chi2","run",condition);
445 SaveGraph(tree,"dcaz_posC_0_Err","run",condition);
446 SaveGraph(tree,"dcaz_posC_1_Err","run",condition);
447 SaveGraph(tree,"dcaz_posC_2_Err","run",condition);
449 SaveGraph(tree,"dcar_posC_0","run",condition);
450 SaveGraph(tree,"dcar_posC_1","run",condition);
451 SaveGraph(tree,"dcar_posC_2","run",condition);
452 SaveGraph(tree,"dcar_posC_chi2","run",condition);
453 SaveGraph(tree,"dcar_posC_0_Err","run",condition);
454 SaveGraph(tree,"dcar_posC_1_Err","run",condition);
455 SaveGraph(tree,"dcar_posC_2_Err","run",condition);
457 SaveGraph(tree,"dcar_negA_0","run",condition);
458 SaveGraph(tree,"dcar_negA_1","run",condition);
459 SaveGraph(tree,"dcar_negA_2","run",condition);
460 SaveGraph(tree,"dcar_negA_chi2","run",condition);
461 SaveGraph(tree,"dcar_negA_0_Err","run",condition);
462 SaveGraph(tree,"dcar_negA_1_Err","run",condition);
463 SaveGraph(tree,"dcar_negA_2_Err","run",condition);
465 SaveGraph(tree,"dcaz_negA_0","run",condition);
466 SaveGraph(tree,"dcaz_negA_1","run",condition);
467 SaveGraph(tree,"dcaz_negA_2","run",condition);
468 SaveGraph(tree,"dcaz_negA_chi2","run",condition);
469 SaveGraph(tree,"dcaz_negA_0_Err","run",condition);
470 SaveGraph(tree,"dcaz_negA_1_Err","run",condition);
471 SaveGraph(tree,"dcaz_negA_2_Err","run",condition);
473 SaveGraph(tree,"dcaz_negC_0","run",condition);
474 SaveGraph(tree,"dcaz_negC_1","run",condition);
475 SaveGraph(tree,"dcaz_negC_2","run",condition);
476 SaveGraph(tree,"dcaz_negC_chi2","run",condition);
477 SaveGraph(tree,"dcaz_negC_0_Err","run",condition);
478 SaveGraph(tree,"dcaz_negC_1_Err","run",condition);
479 SaveGraph(tree,"dcaz_negC_2_Err","run",condition);
481 SaveGraph(tree,"dcar_negC_0","run",condition);
482 SaveGraph(tree,"dcar_negC_1","run",condition);
483 SaveGraph(tree,"dcar_negC_2","run",condition);
484 SaveGraph(tree,"dcar_negC_chi2","run",condition);
485 SaveGraph(tree,"dcar_negC_0_Err","run",condition);
486 SaveGraph(tree,"dcar_negC_1_Err","run",condition);
487 SaveGraph(tree,"dcar_negC_2_Err","run",condition);
491 ////////////////////////////////////////////////////////////////////////////////////////////////////////
497 if (tree) { delete tree; tree=0; }
498 if (out) { delete out; out=0; }
502 //_____________________________________________________________________________
503 Int_t AliTPCPerformanceSummary::SaveGraph(TTree* tree, const Char_t* y, const Char_t* x, const Char_t* condition)
506 // Creates a Graph and writes the canvas to the current directory
507 // called by ProduceTrends function.
513 tree->Draw(s.Data(),condition,"goff");
514 TCanvas* c = new TCanvas(s.Data(),s.Data());
516 TPad* p = new TPad("pad0","pad0",0,0,1,1);
519 if (tree->GetSelectedRows() > 0) {
520 TGraph* graph = new TGraph(tree->GetSelectedRows(), tree->GetV2(), tree->GetV1());
522 graph->GetXaxis()->SetTitle(x);
523 graph->GetYaxis()->SetTitle(y);
527 //graph->Write(s.Data());
533 //_____________________________________________________________________________
534 Int_t AliTPCPerformanceSummary::AnalyzeDCARPhi(const AliPerformanceTPC* pTPC, TTreeSRedirector* const pcstream)
537 // Analyse DCA R imperfections
540 if (!pcstream) return 8;
544 static Double_t offsetdRA=0;
545 static Double_t slopedRA=0;
546 static Double_t offsetdRC=0;
547 static Double_t slopedRC=0;
548 static Double_t offsetdRAErr=0;
549 static Double_t slopedRAErr=0;
550 static Double_t offsetdRCErr=0;
551 static Double_t slopedRCErr=0;
552 static Double_t offsetdRAchi2=0;
553 static Double_t slopedRAchi2=0;
554 static Double_t offsetdRCchi2=0;
555 static Double_t slopedRCchi2=0;
556 static Double_t dcarAP0 = 0;
557 static Double_t dcarAP1 = 0;
558 static Double_t dcarCP0 = 0;
559 static Double_t dcarCP1 = 0;
561 //AliPerformanceTPC* pTPC = dynamic_cast<AliPerformanceTPC*>(pTPCObject);
567 if (pTPC->GetHistos()->FindObject("h_tpc_track_all_recvertex_3_5_7")) {
568 his3D = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_all_recvertex_3_5_7"));
570 his3D->GetYaxis()->SetRangeUser(-1,1);
571 his3D->GetZaxis()->SetRangeUser(0.25,10);
574 static TF1 *fpol1 = new TF1("fpol1","pol1");
576 if (his3D && !fgForceTHnSparse) {
577 his2D = dynamic_cast<TH2*>(his3D->Project3D("xy"));
579 his2D = pTPC->GetTPCTrackHisto()->Projection(3,5);
584 his2D->FitSlicesY(0,0,-1,10,"QNR",&arrayFit);
586 his1D = (TH1*) arrayFit.At(1);
587 his1D->Fit(fpol1,"QNRROB=0.8","QNR",-0.8,-0.1);
588 offsetdRC=fpol1->GetParameter(0);
589 slopedRC=fpol1->GetParameter(1);
590 offsetdRCchi2=fpol1->GetChisquare();
591 slopedRCchi2=fpol1->GetChisquare();
593 his1D->Fit(fpol1,"QNRROB=0.8","QNR",0.1,0.8);
594 offsetdRA=fpol1->GetParameter(0);
595 slopedRA=fpol1->GetParameter(1);
596 offsetdRAErr=fpol1->GetParError(0);
597 slopedRAErr=fpol1->GetParError(1);
598 offsetdRAchi2=fpol1->GetChisquare();
599 slopedRAchi2=fpol1->GetChisquare();
601 printf("DCA R QA report\n");
602 printf("offsetdRA\t%f\n",offsetdRA);
603 printf("slopedRA\t%f\n",slopedRA);
604 printf("offsetdRC\t%f\n",offsetdRC);
605 printf("slopedRC\t%f\n",slopedRC);
608 //extraction of DCAr versus pt
610 TLinearFitter linearFit;
611 linearFit.SetFormula("pol1");
612 TObjArray arrayWidth;
620 his3D->GetYaxis()->SetRangeUser(-1,1);
622 //get his2D in A Side
623 his3D->GetYaxis()->SetRangeUser(0,1);
624 his3D->GetZaxis()->SetRangeUser(0.35,8);
625 his2D = dynamic_cast<TH2*>(his3D->Project3D("xz"));
627 his2D->FitSlicesY(0,0,-1,0,"QNR",&arrayWidth);
629 width = dynamic_cast<TH1*>(arrayWidth.At(2));
632 nXbins = width->GetNbinsX();
633 for(Int_t i=2; i<nXbins; i++){
634 x = width->GetBinCenter(i);
637 y = width->GetBinContent(i);
639 linearFit.AddPoint(&x,y,1);
641 if(!linearFit.Eval()){
643 dcarAP0 = linearFit.GetParameter(0);
645 pn = Int_t(TMath::Abs(dcarAP0)/dcarAP0);
646 dcarAP0 = pn*TMath::Sqrt(TMath::Abs(dcarAP0));
648 dcarAP1 = linearFit.GetParameter(1);
650 pn = Int_t(TMath::Abs(dcarAP1)/dcarAP1);
651 dcarAP1 = pn*TMath::Sqrt(TMath::Abs(dcarAP1));
655 linearFit.ClearPoints();
657 //get his2D in C Side
658 his3D->GetYaxis()->SetRangeUser(-1,-0.001);
659 his2D = dynamic_cast<TH2*>(his3D->Project3D("xz"));
661 his2D->FitSlicesY(0,0,-1,0,"QNR",&arrayWidth);
662 width = dynamic_cast<TH1*>(arrayWidth.At(2));
665 nXbins = width->GetNbinsX();
666 for(Int_t i=2; i<nXbins; i++){
667 x = width->GetBinCenter(i);
670 y = width->GetBinContent(i);
672 linearFit.AddPoint(&x,y);
674 if(!linearFit.Eval()){
675 dcarCP0 = linearFit.GetParameter(0);
677 pn = Int_t(TMath::Abs(dcarCP0)/dcarCP0);
678 dcarCP0 = pn*TMath::Sqrt(TMath::Abs(dcarCP0));
680 dcarCP1 = linearFit.GetParameter(1);
682 pn = Int_t(TMath::Abs(dcarCP1)/dcarCP1);
683 dcarCP1 = pn*TMath::Sqrt(TMath::Abs(dcarCP1));
686 his3D->GetYaxis()->SetRangeUser(-1,1);
687 his3D->GetZaxis()->SetRangeUser(0,20);
692 (*pcstream)<<"tpcQA"<<
693 "offsetdRA="<< offsetdRA<< // mean r-phi distortion fit A side (DCA_rphi=[0]+[1]*tgl(theta)) - parameter [0] offset
694 "slopedRA="<< slopedRA<< // mean r-phi distortion fit A side (DCA_rphi=[0]+[1]*tgl(theta)) - parameter [1] slope
695 "offsetdRC="<< offsetdRC<< //
696 "slopedRC="<<slopedRC<<
698 "offsetdRAErr="<< offsetdRAErr<<
699 "slopedRAErr="<< slopedRAErr<<
700 "offsetdRCErr="<< offsetdRCErr<<
701 "slopedRCErr="<<slopedRCErr<<
703 "offsetdRAchi2="<< offsetdRAchi2<<
704 "slopedRAchi2="<< slopedRAchi2<<
705 "offsetdRCchi2="<< offsetdRCchi2<<
706 "slopedRCchi2="<<slopedRCchi2<<
708 "dcarAP0="<<dcarAP0<<
709 "dcarAP1="<<dcarAP1<<
710 "dcarCP0="<<dcarCP0<<
716 //_____________________________________________________________________________
717 Int_t AliTPCPerformanceSummary::AnalyzeDCARPhiPos(const AliPerformanceTPC* pTPC, TTreeSRedirector* const pcstream)
720 // Analyse DCA R imperfections for positive particles
723 if (!pcstream) return 16;
724 if (!pTPC) return 16;
727 static Double_t offsetdRAPos=0;
728 static Double_t slopedRAPos=0;
729 static Double_t offsetdRCPos=0;
730 static Double_t slopedRCPos=0;
731 static Double_t offsetdRAErrPos=0;
732 static Double_t slopedRAErrPos=0;
733 static Double_t offsetdRCErrPos=0;
734 static Double_t slopedRCErrPos=0;
735 static Double_t offsetdRAchi2Pos=0;
736 static Double_t slopedRAchi2Pos=0;
737 static Double_t offsetdRCchi2Pos=0;
738 static Double_t slopedRCchi2Pos=0;
740 //AliPerformanceTPC* pTPC = dynamic_cast<AliPerformanceTPC*>(pTPCObject);
745 if (pTPC->GetHistos()->FindObject("h_tpc_track_pos_recvertex_3_5_7")) {
746 his3D = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_pos_recvertex_3_5_7"));
747 if(!his3D) return 16;
748 his3D->GetYaxis()->SetRangeUser(-1,1);
749 his3D->GetZaxis()->SetRangeUser(0.25,10);
752 static TF1 *fpol1 = new TF1("fpol1","pol1");
754 if (his3D && !fgForceTHnSparse) {
755 his2D = dynamic_cast<TH2*>(his3D->Project3D("xy"));
757 pTPC->GetTPCTrackHisto()->GetAxis(8)->SetRangeUser(0,1.5);
758 his2D = pTPC->GetTPCTrackHisto()->Projection(3,5);
759 pTPC->GetTPCTrackHisto()->GetAxis(8)->SetRangeUser(-1.5,1.5);
761 if(!his2D) return 16;
763 his2D->FitSlicesY(0,0,-1,10,"QNR",&arrayFit);
765 his1D = (TH1*) arrayFit.At(1);
766 his1D->Fit(fpol1,"QNRROB=0.8","QNR",-0.8,-0.1);
767 offsetdRCPos=fpol1->GetParameter(0);
768 slopedRCPos=fpol1->GetParameter(1);
769 offsetdRCchi2Pos=fpol1->GetChisquare();
770 slopedRCchi2Pos=fpol1->GetChisquare();
772 his1D->Fit(fpol1,"QNRROB=0.8","QNR",0.1,0.8);
773 offsetdRAPos=fpol1->GetParameter(0);
774 slopedRAPos=fpol1->GetParameter(1);
775 offsetdRAErrPos=fpol1->GetParError(0);
776 slopedRAErrPos=fpol1->GetParError(1);
777 offsetdRAchi2Pos=fpol1->GetChisquare();
778 slopedRAchi2Pos=fpol1->GetChisquare();
780 printf("DCA R QA Pos report\n");
781 printf("offsetdRAPos\t%f\n",offsetdRAPos);
782 printf("slopedRAPos\t%f\n",slopedRAPos);
783 printf("offsetdRCPos\t%f\n",offsetdRCPos);
784 printf("slopedRCPos\t%f\n",slopedRCPos);
788 (*pcstream)<<"tpcQA"<<
789 "offsetdRAPos="<< offsetdRAPos<<
790 "slopedRAPos="<< slopedRAPos<<
791 "offsetdRCPos="<< offsetdRCPos<<
792 "slopedRCPos="<<slopedRCPos<<
794 "offsetdRAErrPos="<< offsetdRAErrPos<<
795 "slopedRAErrPos="<< slopedRAErrPos<<
796 "offsetdRCErrPos="<< offsetdRCErrPos<<
797 "slopedRCErrPos="<<slopedRCErrPos<<
799 "offsetdRAchi2Pos="<< offsetdRAchi2Pos<<
800 "slopedRAchi2Pos="<< slopedRAchi2Pos<<
801 "offsetdRCchi2Pos="<< offsetdRCchi2Pos<<
802 "slopedRCchi2Pos="<<slopedRCchi2Pos;
807 //_____________________________________________________________________________
808 Int_t AliTPCPerformanceSummary::AnalyzeDCARPhiNeg(const AliPerformanceTPC* pTPC, TTreeSRedirector* const pcstream)
811 // Analyse DCA R imperfections for negative particles
813 if (!pcstream) return 32;
814 if (!pTPC) return 32;
817 static Double_t offsetdRANeg=0;
818 static Double_t slopedRANeg=0;
819 static Double_t offsetdRCNeg=0;
820 static Double_t slopedRCNeg=0;
821 static Double_t offsetdRAErrNeg=0;
822 static Double_t slopedRAErrNeg=0;
823 static Double_t offsetdRCErrNeg=0;
824 static Double_t slopedRCErrNeg=0;
825 static Double_t offsetdRAchi2Neg=0;
826 static Double_t slopedRAchi2Neg=0;
827 static Double_t offsetdRCchi2Neg=0;
828 static Double_t slopedRCchi2Neg=0;
830 //AliPerformanceTPC* pTPC = dynamic_cast<AliPerformanceTPC*>(pTPCObject);
835 if (pTPC->GetHistos()->FindObject("h_tpc_track_neg_recvertex_3_5_7")) {
836 his3D = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_neg_recvertex_3_5_7"));
837 if(!his3D) return 32;
838 his3D->GetYaxis()->SetRangeUser(-1,1);
839 his3D->GetZaxis()->SetRangeUser(0.25,10);
842 static TF1 *fpol1 = new TF1("fpol1","pol1");
844 if (his3D && !fgForceTHnSparse) {
845 his2D = dynamic_cast<TH2*>(his3D->Project3D("xy"));
847 pTPC->GetTPCTrackHisto()->GetAxis(8)->SetRangeUser(-1.5,0);
848 his2D = pTPC->GetTPCTrackHisto()->Projection(3,5);
849 pTPC->GetTPCTrackHisto()->GetAxis(8)->SetRangeUser(-1.5,1.5);
851 if(!his2D) return 32;
853 his2D->FitSlicesY(0,0,-1,10,"QNR",&arrayFit);
855 his1D = (TH1*) arrayFit.At(1);
856 his1D->Fit(fpol1,"QNRROB=0.8","QNR",-0.8,-0.1);
857 offsetdRCNeg=fpol1->GetParameter(0);
858 slopedRCNeg=fpol1->GetParameter(1);
859 offsetdRCchi2Neg=fpol1->GetChisquare();
860 slopedRCchi2Neg=fpol1->GetChisquare();
862 his1D->Fit(fpol1,"QNRROB=0.8","QNR",0.1,0.8);
863 offsetdRANeg=fpol1->GetParameter(0);
864 slopedRANeg=fpol1->GetParameter(1);
865 offsetdRAErrNeg=fpol1->GetParError(0);
866 slopedRAErrNeg=fpol1->GetParError(1);
867 offsetdRAchi2Neg=fpol1->GetChisquare();
868 slopedRAchi2Neg=fpol1->GetChisquare();
870 printf("DCA R QA Neg report\n");
871 printf("offsetdRANeg\t%f\n",offsetdRANeg);
872 printf("slopedRANeg\t%f\n",slopedRANeg);
873 printf("offsetdRCNeg\t%f\n",offsetdRCNeg);
874 printf("slopedRCNeg\t%f\n",slopedRCNeg);
876 // dump drift QA values
878 (*pcstream)<<"tpcQA"<<
879 "offsetdRANeg="<< offsetdRANeg<<
880 "slopedRANeg="<< slopedRANeg<<
881 "offsetdRCNeg="<< offsetdRCNeg<<
882 "slopedRCNeg="<<slopedRCNeg<<
884 "offsetdRAErrNeg="<< offsetdRAErrNeg<<
885 "slopedRAErrNeg="<< slopedRAErrNeg<<
886 "offsetdRCErrNeg="<< offsetdRCErrNeg<<
887 "slopedRCErrNeg="<<slopedRCErrNeg<<
889 "offsetdRAchi2Neg="<< offsetdRAchi2Neg<<
890 "slopedRAchi2Neg="<< slopedRAchi2Neg<<
891 "offsetdRCchi2Neg="<< offsetdRCchi2Neg<<
892 "slopedRCchi2Neg="<<slopedRCchi2Neg;
897 //_____________________________________________________________________________
898 Int_t AliTPCPerformanceSummary::AnalyzeNCL(const AliPerformanceTPC* pTPC, TTreeSRedirector* const pcstream)
901 // Analyse number of TPC clusters
904 if (!pcstream) return 1;
908 static Double_t meanTPCnclF=0;
909 static Double_t meanTPCnclFStat=0;
910 static Double_t meanTPCnclFErr=0;
911 static Double_t rmsTPCnclF=0;
912 static Double_t meanTPCChi2=0;
913 static Double_t rmsTPCChi2=0;
914 static Double_t slopeATPCnclF=0;
915 static Double_t slopeCTPCnclF=0;
916 static Double_t slopeATPCnclFErr=0;
917 static Double_t slopeCTPCnclFErr=0;
918 static Double_t meanTPCncl=0;
919 static Double_t meanTPCnclStat=0;
920 static Double_t meanTPCnclErr=0;
921 static Double_t rmsTPCncl=0;
922 static Double_t slopeATPCncl=0;
923 static Double_t slopeCTPCncl=0;
924 static Double_t slopeATPCnclErr=0;
925 static Double_t slopeCTPCnclErr=0;
932 static TF1 *fpol1 = new TF1("fpol1","pol1");
935 // only events with rec. vertex
937 // pt cut - 0.250 GeV
938 pTPC->GetTPCTrackHisto()->GetAxis(5)->SetRangeUser(-1.,1.);
939 pTPC->GetTPCTrackHisto()->GetAxis(7)->SetRangeUser(0.25,10);
941 if (pTPC->GetHistos()->FindObject("h_tpc_track_all_recvertex_0_5_7")) {
942 his3D0 = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_all_recvertex_0_5_7"));
943 if(!his3D0) return 1;
944 his3D0->GetYaxis()->SetRangeUser(-1,1);
945 his3D0->GetZaxis()->SetRangeUser(0.25,10);
947 if (pTPC->GetHistos()->FindObject("h_tpc_track_all_recvertex_1_5_7")) {
948 his3D1 = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_all_recvertex_1_5_7"));
949 if(!his3D1) return 1;
950 his3D1->GetYaxis()->SetRangeUser(-1,1);
951 his3D1->GetZaxis()->SetRangeUser(0.25,10);
953 if (pTPC->GetHistos()->FindObject("h_tpc_track_all_recvertex_2_5_7")) {
954 his3D2 = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_all_recvertex_2_5_7"));
955 if(!his3D2) return 1;
956 his3D2->GetYaxis()->SetRangeUser(-1,1);
957 his3D2->GetZaxis()->SetRangeUser(0.25,10);
958 his3D2->GetXaxis()->SetRangeUser(0.4,1.1);
962 if (his3D0 && !fgForceTHnSparse) {
963 his1D = his3D0->Project3D("x");
965 his1D = pTPC->GetTPCTrackHisto()->Projection(0);
968 meanTPCncl= his1D->GetMean();
969 meanTPCnclStat= his1D->GetEntries();
970 meanTPCnclErr= his1D->GetMeanError();
971 rmsTPCncl= his1D->GetRMS();
975 if (his3D1 && !fgForceTHnSparse) {
976 his1D = his3D1->Project3D("x");
978 his1D = pTPC->GetTPCTrackHisto()->Projection(1);
981 meanTPCChi2= his1D->GetMean();
982 rmsTPCChi2= his1D->GetRMS();
985 if (his3D0 && !fgForceTHnSparse) {
986 hprof = (dynamic_cast<TH2*>(his3D0->Project3D("xy")))->ProfileX();
988 hprof = pTPC->GetTPCTrackHisto()->Projection(0,5)->ProfileX();
992 hprof->Fit(fpol1,"QNR","QNR",0.1,0.8);
993 slopeATPCncl= fpol1->GetParameter(1);
994 slopeATPCnclErr= fpol1->GetParError(1);
995 hprof->Fit(fpol1,"QNR","QNR",-0.8,-0.1);
996 slopeCTPCncl= fpol1->GetParameter(1);
997 slopeCTPCnclErr= fpol1->GetParameter(1);
1001 // findable clusters
1004 if (his3D2 && !fgForceTHnSparse) {
1005 his1D = his3D2->Project3D("x");
1007 pTPC->GetTPCTrackHisto()->GetAxis(2)->SetRangeUser(0.4,1.1);
1008 his1D = pTPC->GetTPCTrackHisto()->Projection(2);
1011 meanTPCnclF= his1D->GetMean();
1012 meanTPCnclFStat= his1D->GetEntries();
1013 meanTPCnclFErr= his1D->GetMeanError();
1014 rmsTPCnclF= his1D->GetRMS();
1017 if (his3D2 && !fgForceTHnSparse) {
1018 his1D = (dynamic_cast<TH2*>(his3D2->Project3D("xy")))->ProfileX();
1020 pTPC->GetTPCTrackHisto()->GetAxis(2)->SetRangeUser(0.4,1.1);
1021 his1D = pTPC->GetTPCTrackHisto()->Projection(2,5)->ProfileX();
1023 if(!his1D) return 1;
1025 his1D->Fit(fpol1,"QNR","QNR",0.1,0.8);
1026 slopeATPCnclF= fpol1->GetParameter(1);
1027 slopeATPCnclFErr= fpol1->GetParError(1);
1028 his1D->Fit(fpol1,"QNR","QNR",-0.8,-0.1);
1029 slopeCTPCnclF= fpol1->GetParameter(1);
1030 slopeCTPCnclFErr= fpol1->GetParameter(1);
1033 pTPC->GetTPCTrackHisto()->GetAxis(2)->SetRangeUser(0,10);
1035 printf("Cluster QA report\n");
1036 printf("meanTPCnclF=\t%f\n",meanTPCnclF);
1037 printf("rmsTPCnclF=\t%f\n",rmsTPCnclF);
1038 printf("slopeATPCnclF=\t%f\n",slopeATPCnclF);
1039 printf("slopeCTPCnclF=\t%f\n",slopeCTPCnclF);
1040 printf("meanTPCncl=\t%f\n",meanTPCncl);
1041 printf("rmsTPCncl=\t%f\n",rmsTPCncl);
1042 printf("slopeATPCncl=\t%f\n",slopeATPCncl);
1043 printf("slopeCTPCncl=\t%f\n",slopeCTPCncl);
1044 printf("meanTPCChi2=\t%f\n",meanTPCChi2);
1045 printf("rmsTPCChi2=\t%f\n",rmsTPCChi2);
1047 // dump results to the tree
1049 (*pcstream)<<"tpcQA"<<
1050 "meanTPCnclF="<<meanTPCnclF << // mean number of cluster/number of findable cluster
1051 "meanTPCnclFStat="<<meanTPCnclFStat << // number of cluster/number of findable cluster - number of etries
1052 "meanTPCnclFErr="<<meanTPCnclFErr << // error of mean number of cluster/number of findable cluster
1053 "rmsTPCnclF="<<rmsTPCnclF << // RMS of distribution of number of cluster/number of findable cluster
1054 "meanTPCChi2="<<meanTPCChi2 << // ????
1055 "rmsTPCChi2="<<rmsTPCChi2 << // ????
1056 "slopeATPCnclF="<< slopeATPCnclF<< // A side- fit of number of clusters/findable (Ncl=[0]+[1]*tan(\theta)) - parameter [1] slope
1057 "slopeCTPCnclF="<< slopeCTPCnclF<< // C side- fit of number of clusters/findable (Ncl=[0]+[1]*tan(\theta)) - parameter [1] slope
1058 "slopeATPCnclFErr="<< slopeATPCnclFErr<< // A side- fit of number of clusters/findable (Ncl=[0]+[1]*tan(\theta)) - error of parameter [1] slope
1059 "slopeCTPCnclFErr="<< slopeCTPCnclFErr<< // C side- fit of number of clusters/findable (Ncl=[0]+[1]*tan(\theta)) - error of parameter [1] slope
1060 "meanTPCncl="<<meanTPCncl << // mean number of cluster
1061 "meanTPCnclStat="<<meanTPCnclStat << // number of cluster - number of etries in histogram
1062 "meanTPCnclErr="<<meanTPCnclErr << // error of mean number of cluster
1063 "rmsTPCncl="<< rmsTPCncl<< // rms of mean number of cluster
1064 "slopeATPCncl="<< slopeATPCncl<< // A side- fit of number of clusters (Ncl=[0]+[1]*tan(\theta)) - parameter [1] slope
1065 "slopeCTPCncl="<< slopeCTPCncl<< // C side - fit of number of clusters (Ncl=[0]+[1]*tan(\theta)) - parameter 1 - slope
1066 "slopeATPCnclErr="<< slopeATPCnclErr<< // A side - fit of number of clusters (Ncl=[0]+[1]*tan(\theta)) - error of parameter 1 - slope
1067 "slopeCTPCnclErr="<< slopeCTPCnclErr; // C side - fit of number of clusters (Ncl=[0]+[1]*tan(\theta)) - error of parameter 1 - slope
1072 //_____________________________________________________________________________
1073 Int_t AliTPCPerformanceSummary::AnalyzeDrift(const AliPerformanceTPC* pTPC, TTreeSRedirector* const pcstream)
1076 // Analyse DCA Z imperferctions (drift velocity)
1079 if (!pcstream) return 2;
1080 if (!pTPC) return 2;
1083 static Double_t offsetdZA=0;
1084 static Double_t slopedZA=0;
1085 static Double_t offsetdZC=0;
1086 static Double_t slopedZC=0;
1087 static Double_t offsetdZAErr=0;
1088 static Double_t slopedZAErr=0;
1089 static Double_t offsetdZCErr=0;
1090 static Double_t slopedZCErr=0;
1091 static Double_t offsetdZAchi2=0;
1092 static Double_t slopedZAchi2=0;
1093 static Double_t offsetdZCchi2=0;
1094 static Double_t slopedZCchi2=0;
1099 if (pTPC->GetHistos()->FindObject("h_tpc_track_all_recvertex_4_5_7")) {
1100 his3D = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_all_recvertex_4_5_7"));
1101 if(!his3D) return 2;
1102 his3D->GetYaxis()->SetRangeUser(-1,1);
1103 his3D->GetZaxis()->SetRangeUser(0.25,10);
1106 if (his3D && !fgForceTHnSparse) {
1107 his2D = dynamic_cast<TH2*>(his3D->Project3D("xy"));
1109 his2D = pTPC->GetTPCTrackHisto()->Projection(4,5);
1111 if(!his2D) return 2;
1113 static TF1 *fpol1 = new TF1("fpol1","pol1");
1115 his2D->FitSlicesY(0,0,-1,10,"QNR",&arrayFit);
1117 his1D = (TH1*) arrayFit.At(1);
1118 his1D->Fit(fpol1,"QNRROB=0.8","QNR",-0.8,-0.1);
1119 offsetdZC=fpol1->GetParameter(0);
1120 slopedZC=fpol1->GetParameter(1);
1121 offsetdZCErr=fpol1->GetParError(0);
1122 slopedZCErr=fpol1->GetParError(1);
1123 offsetdZCchi2=fpol1->GetChisquare();
1124 slopedZCchi2=fpol1->GetChisquare();
1126 his1D->Fit(fpol1,"QNRROB=0.8","QNR",0.1,0.8);
1127 offsetdZA=fpol1->GetParameter(0);
1128 slopedZA=fpol1->GetParameter(1);
1129 offsetdZAErr=fpol1->GetParError(0);
1130 slopedZAErr=fpol1->GetParError(1);
1131 offsetdZAchi2=fpol1->GetChisquare();
1132 slopedZAchi2=fpol1->GetChisquare();
1134 printf("Drift velocity QA report\n");
1135 printf("offsetdZA\t%f\n",offsetdZA);
1136 printf("slopedZA\t%f\n",slopedZA);
1137 printf("offsetdZC\t%f\n",offsetdZC);
1138 printf("slopedZC\t%f\n",slopedZC);
1140 // dump drift QA values
1142 (*pcstream)<<"tpcQA"<<
1143 "offsetdZA="<< offsetdZA<< //
1144 "slopedZA="<< slopedZA<<
1145 "offsetdZC="<< offsetdZC<<
1146 "slopedZC="<<slopedZC<<
1148 "offsetdZAErr="<< offsetdZAErr<<
1149 "slopedZAErr="<< slopedZAErr<<
1150 "offsetdZCErr="<< offsetdZCErr<<
1151 "slopedZCErr="<<slopedZCErr<<
1153 "offsetdZAchi2="<< offsetdZAchi2<<
1154 "slopedZAchi2="<< slopedZAchi2<<
1155 "offsetdZCchi2="<< offsetdZCchi2<<
1156 "slopedZCchi2="<<slopedZCchi2;
1161 //_____________________________________________________________________________
1162 Int_t AliTPCPerformanceSummary::AnalyzeDriftPos(const AliPerformanceTPC* pTPC, TTreeSRedirector* const pcstream)
1165 // Analyse DCA Z imperferctions (drift velocity)
1166 // for positive particles
1168 if (!pcstream) return 64;
1169 if (!pTPC) return 64;
1172 static Double_t offsetdZAPos=0;
1173 static Double_t slopedZAPos=0;
1174 static Double_t offsetdZCPos=0;
1175 static Double_t slopedZCPos=0;
1176 static Double_t offsetdZAErrPos=0;
1177 static Double_t slopedZAErrPos=0;
1178 static Double_t offsetdZCErrPos=0;
1179 static Double_t slopedZCErrPos=0;
1180 static Double_t offsetdZAchi2Pos=0;
1181 static Double_t slopedZAchi2Pos=0;
1182 static Double_t offsetdZCchi2Pos=0;
1183 static Double_t slopedZCchi2Pos=0;
1188 if (pTPC->GetHistos()->FindObject("h_tpc_track_pos_recvertex_4_5_7")) {
1189 his3D = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_pos_recvertex_4_5_7"));
1190 if(!his3D) return 64;
1191 his3D->GetYaxis()->SetRangeUser(-1,1);
1192 his3D->GetZaxis()->SetRangeUser(0.25,10);
1195 if (his3D && !fgForceTHnSparse) {
1196 his2D = dynamic_cast<TH2*>(his3D->Project3D("xy"));
1198 his2D = pTPC->GetTPCTrackHisto()->Projection(4,5);
1200 if(!his2D) return 64;
1202 static TF1 *fpol1 = new TF1("fpol1","pol1");
1204 his2D->FitSlicesY(0,0,-1,10,"QNR",&arrayFit);
1207 his1D = (TH1*) arrayFit.At(1);
1208 his1D->Fit(fpol1,"QNRROB=0.8","QNR",-0.8,-0.1);
1209 offsetdZCPos=fpol1->GetParameter(0);
1210 slopedZCPos=fpol1->GetParameter(1);
1211 offsetdZCErrPos=fpol1->GetParError(0);
1212 slopedZCErrPos=fpol1->GetParError(1);
1213 offsetdZCchi2Pos=fpol1->GetChisquare();
1214 slopedZCchi2Pos=fpol1->GetChisquare();
1216 his1D->Fit(fpol1,"QNRROB=0.8","QNR",0.1,0.8);
1217 offsetdZAPos=fpol1->GetParameter(0);
1218 slopedZAPos=fpol1->GetParameter(1);
1219 offsetdZAErrPos=fpol1->GetParError(0);
1220 slopedZAErrPos=fpol1->GetParError(1);
1221 offsetdZAchi2Pos=fpol1->GetChisquare();
1222 slopedZAchi2Pos=fpol1->GetChisquare();
1224 printf("Drift velocity QA report\n");
1225 printf("offsetdZAPos\t%f\n",offsetdZAPos);
1226 printf("slopedZAPos\t%f\n",slopedZAPos);
1227 printf("offsetdZCPos\t%f\n",offsetdZCPos);
1228 printf("slopedZCPos\t%f\n",slopedZCPos);
1230 // dump drift QA values
1232 (*pcstream)<<"tpcQA"<<
1233 "offsetdZAPos="<< offsetdZAPos<<
1234 "slopedZAPos="<< slopedZAPos<<
1235 "offsetdZCPos="<< offsetdZCPos<<
1236 "slopedZCPos="<<slopedZCPos<<
1238 "offsetdZAErrPos="<< offsetdZAErrPos<<
1239 "slopedZAErrPos="<< slopedZAErrPos<<
1240 "offsetdZCErrPos="<< offsetdZCErrPos<<
1241 "slopedZCErrPos="<<slopedZCErrPos<<
1243 "offsetdZAchi2Pos="<< offsetdZAchi2Pos<<
1244 "slopedZAchi2Pos="<< slopedZAchi2Pos<<
1245 "offsetdZCchi2Pos="<< offsetdZCchi2Pos<<
1246 "slopedZCchi2Pos="<<slopedZCchi2Pos;
1251 //_____________________________________________________________________________
1252 Int_t AliTPCPerformanceSummary::AnalyzeDriftNeg(const AliPerformanceTPC* pTPC, TTreeSRedirector* const pcstream)
1255 // Analyse DCA Z imperferctions (drift velocity)
1256 // for negative particles
1258 if (!pcstream) return 128;
1259 if (!pTPC) return 128;
1262 static Double_t offsetdZANeg=0;
1263 static Double_t slopedZANeg=0;
1264 static Double_t offsetdZCNeg=0;
1265 static Double_t slopedZCNeg=0;
1266 static Double_t offsetdZAErrNeg=0;
1267 static Double_t slopedZAErrNeg=0;
1268 static Double_t offsetdZCErrNeg=0;
1269 static Double_t slopedZCErrNeg=0;
1270 static Double_t offsetdZAchi2Neg=0;
1271 static Double_t slopedZAchi2Neg=0;
1272 static Double_t offsetdZCchi2Neg=0;
1273 static Double_t slopedZCchi2Neg=0;
1279 if (pTPC->GetHistos()->FindObject("h_tpc_track_neg_recvertex_4_5_7")) {
1280 his3D = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_neg_recvertex_4_5_7"));
1281 if(!his3D) return 128;
1282 his3D->GetYaxis()->SetRangeUser(-1,1);
1283 his3D->GetZaxis()->SetRangeUser(0.25,10);
1285 if (his3D && !fgForceTHnSparse) {
1286 his2D = dynamic_cast<TH2*>(his3D->Project3D("xy"));
1288 his2D = pTPC->GetTPCTrackHisto()->Projection(4,5);
1290 if(!his2D) return 128;
1292 static TF1 *fpol1 = new TF1("fpol1","pol1");
1294 his2D->FitSlicesY(0,0,-1,10,"QNR",&arrayFit);
1297 his1D = (TH1*) arrayFit.At(1);
1298 his1D->Fit(fpol1,"QNRROB=0.8","QNR",-0.8,-0.1);
1299 offsetdZCNeg=fpol1->GetParameter(0);
1300 slopedZCNeg=fpol1->GetParameter(1);
1301 offsetdZCErrNeg=fpol1->GetParError(0);
1302 slopedZCErrNeg=fpol1->GetParError(1);
1303 offsetdZCchi2Neg=fpol1->GetChisquare();
1304 slopedZCchi2Neg=fpol1->GetChisquare();
1306 his1D->Fit(fpol1,"QNRROB=0.8","QNR",0.1,0.8);
1307 offsetdZANeg=fpol1->GetParameter(0);
1308 slopedZANeg=fpol1->GetParameter(1);
1309 offsetdZAErrNeg=fpol1->GetParError(0);
1310 slopedZAErrNeg=fpol1->GetParError(1);
1311 offsetdZAchi2Neg=fpol1->GetChisquare();
1312 slopedZAchi2Neg=fpol1->GetChisquare();
1314 printf("Drift velocity QA report\n");
1315 printf("offsetdZANeg\t%f\n",offsetdZANeg);
1316 printf("slopedZANeg\t%f\n",slopedZANeg);
1317 printf("offsetdZCNeg\t%f\n",offsetdZCNeg);
1318 printf("slopedZCNeg\t%f\n",slopedZCNeg);
1320 // dump drift QA values
1322 (*pcstream)<<"tpcQA"<<
1323 "offsetdZANeg="<< offsetdZANeg<<
1324 "slopedZANeg="<< slopedZANeg<<
1325 "offsetdZCNeg="<< offsetdZCNeg<<
1326 "slopedZCNeg="<<slopedZCNeg<<
1328 "offsetdZAErrNeg="<< offsetdZAErrNeg<<
1329 "slopedZAErrNeg="<< slopedZAErrNeg<<
1330 "offsetdZCErrNeg="<< offsetdZCErrNeg<<
1331 "slopedZCErrNeg="<<slopedZCErrNeg<<
1333 "offsetdZAchi2Neg="<< offsetdZAchi2Neg<<
1334 "slopedZAchi2Neg="<< slopedZAchi2Neg<<
1335 "offsetdZCchi2Neg="<< offsetdZCchi2Neg<<
1336 "slopedZCchi2Neg="<<slopedZCchi2Neg;
1341 //_____________________________________________________________________________
1342 Int_t AliTPCPerformanceSummary::AnalyzeGain(const AliPerformanceDEdx* pTPCgain, TTreeSRedirector* const pcstream)
1348 if (!pcstream) return 4;
1349 if (!pTPCgain) return 4;
1351 static TVectorD meanMIPvsSector(36);
1352 static TVectorD sector(36);
1353 static Float_t meanMIP = 0;
1354 static Float_t resolutionMIP = 0;
1355 static Float_t attachSlopeC = 0;
1356 static Float_t attachSlopeA = 0;
1359 //TH1 * hisProj1D=0;
1363 meanMIPvsSector.Zero();
1365 // select MIP particles
1367 pTPCgain->GetDeDxHisto()->GetAxis(7)->SetRangeUser(0.4,0.55);
1368 pTPCgain->GetDeDxHisto()->GetAxis(0)->SetRangeUser(35,60);
1369 pTPCgain->GetDeDxHisto()->GetAxis(6)->SetRangeUser(80,160);
1370 pTPCgain->GetDeDxHisto()->GetAxis(5)->SetRangeUser(-1,1);
1372 // MIP position and resolution
1374 TF1 gausFit("gausFit","gaus");
1376 if (pTPCgain->GetHistos()->FindObject("h_tpc_dedx_mips_0") && !fgForceTHnSparse) {
1377 his1D = dynamic_cast<TH1*>(pTPCgain->GetHistos()->FindObject("h_tpc_dedx_mips_0")->Clone());
1379 his1D = pTPCgain->GetDeDxHisto()->Projection(0);
1381 if(!his1D) return 4;
1382 his1D->Fit(&gausFit,"QN","QN");
1384 meanMIP = gausFit.GetParameter(1);
1386 if (meanMIP!=0) resolutionMIP = gausFit.GetParameter(2)/meanMIP;
1387 //removedtotest// delete his1D;
1389 // MIP position vs. dip angle (attachment)
1391 pTPCgain->GetDeDxHisto()->GetAxis(5)->SetRangeUser(-3,0); // C side
1392 if (pTPCgain->GetHistos()->FindObject("h_tpc_dedx_mips_c_0_5") && !fgForceTHnSparse) {
1393 his2D = dynamic_cast<TH2*>(pTPCgain->GetHistos()->FindObject("h_tpc_dedx_mips_c_0_5")->Clone());
1395 his2D = pTPCgain->GetDeDxHisto()->Projection(0,5);
1397 if(!his2D) return 4;
1399 TF1 * fpol = new TF1("fpol","pol1");
1401 his2D->FitSlicesY(0,0,-1,10,"QN",&arrayFit);
1402 his1D = (TH1*) arrayFit.At(1);
1403 his1D->Fit(fpol,"QNROB=0.8","QNR",-1,0);
1404 attachSlopeC = fpol->GetParameter(1);
1405 //removedtotest// delete his2D;
1406 //removedtotest// delete his1D;
1408 pTPCgain->GetDeDxHisto()->GetAxis(5)->SetRangeUser(0,3); // A side
1409 if (pTPCgain->GetHistos()->FindObject("h_tpc_dedx_mips_a_0_5") && !fgForceTHnSparse) {
1410 his2D = dynamic_cast<TH2*>(pTPCgain->GetHistos()->FindObject("h_tpc_dedx_mips_a_0_5")->Clone());
1412 his2D = pTPCgain->GetDeDxHisto()->Projection(0,5);
1414 if(!his2D) return 4;
1416 TF1 * fpolA = new TF1("fpolA","pol1");
1417 TObjArray arrayFitA;
1418 his2D->FitSlicesY(0,0,-1,10,"QN",&arrayFit);
1419 his1D = (TH1*) arrayFit.At(1);
1420 his1D->Fit(fpolA,"QNROB=0.8","QN",0,1);
1421 attachSlopeA = fpolA->GetParameter(1);
1422 //removedtotest// delete his2D;
1423 //removedtotest// delete his1D;
1425 // MIP position vs. sector
1427 pTPCgain->GetDeDxHisto()->GetAxis(5)->SetRangeUser(-3,0); // C side
1428 if (pTPCgain->GetHistos()->FindObject("h_tpc_dedx_mips_c_0_1") && !fgForceTHnSparse) {
1429 his2D = dynamic_cast<TH2*>(pTPCgain->GetHistos()->FindObject("h_tpc_dedx_mips_c_0_1")->Clone());
1431 his2D = pTPCgain->GetDeDxHisto()->Projection(0,1);
1433 if(!his2D) return 4;
1435 for(Int_t i = 0; i < 18; i++) { // loop over sectors; correct mapping to be checked!
1437 Float_t phiLow = -TMath::Pi() + i*(20./360.)*(2*TMath::Pi());
1438 Float_t phiUp = -TMath::Pi() + (i+1)*(20./360.)*(2*TMath::Pi());
1439 //pTPCgain->GetDeDxHisto()->GetAxis(1)->SetRangeUser(phiLow,phiUp);
1440 his2D->GetXaxis()->SetRangeUser(phiLow,phiUp);
1441 //his1D = pTPCgain->GetDeDxHisto()->Projection(0);
1442 his1D = his2D->ProjectionY();
1443 TF1 gausFunc("gausFunc","gaus");
1444 his1D->Fit(&gausFunc, "QN");
1445 meanMIPvsSector(i) = gausFunc.GetParameter(1);
1447 //removedtotest// delete his1D;
1449 //removedtotest// delete his2D;
1451 pTPCgain->GetDeDxHisto()->GetAxis(5)->SetRangeUser(0,3); // A side
1452 if (pTPCgain->GetHistos()->FindObject("h_tpc_dedx_mips_a_0_1") && !fgForceTHnSparse) {
1453 his2D = dynamic_cast<TH2*>(pTPCgain->GetHistos()->FindObject("h_tpc_dedx_mips_a_0_1")->Clone());
1455 his2D = pTPCgain->GetDeDxHisto()->Projection(0,1);
1457 if(!his2D) return 4;
1459 for(Int_t i = 0; i < 18; i++) { // loop over sectors; correct mapping to be checked!
1461 Float_t phiLow = -TMath::Pi() + i*(20./360.)*(2*TMath::Pi());
1462 Float_t phiUp = -TMath::Pi() + (i+1)*(20./360.)*(2*TMath::Pi());
1463 //pTPCgain->GetDeDxHisto()->GetAxis(1)->SetRangeUser(phiLow,phiUp);
1464 his2D->GetXaxis()->SetRangeUser(phiLow,phiUp);
1465 //his1D = pTPCgain->GetDeDxHisto()->Projection(0);
1466 his1D = his2D->ProjectionY();
1467 TF1 gausFunc("gausFunc","gaus");
1468 his1D->Fit(&gausFunc, "QN");
1469 meanMIPvsSector(i+18) = gausFunc.GetParameter(1);
1471 //removedtotest// delete his1D;
1473 //removedtotest// delete his2D;
1475 printf("Gain QA report\n");
1476 printf("MIP mean\t%f\n",meanMIP);
1477 printf("MIP resolution\t%f\n",resolutionMIP);
1478 printf("MIPslopeA\t%f\n",attachSlopeA);
1479 printf("MIPslopeC\t%f\n",attachSlopeC);
1482 (*pcstream)<<"tpcQA"<<
1483 "MIPattachSlopeC="<<attachSlopeC<<
1484 "MIPattachSlopeA="<<attachSlopeA<<
1485 "resolutionMIP="<<resolutionMIP<<
1486 "meanMIPvsSector.="<<&meanMIPvsSector<<
1487 "sector.="<<§or<<
1488 "meanMIP="<<meanMIP;
1493 //_____________________________________________________________________________
1494 Int_t AliTPCPerformanceSummary::AnalyzeEvent(const AliPerformanceTPC* pTPC, TTreeSRedirector* const pcstream)
1497 // Analyse Primary Vertex Distribution and Multiplicities
1499 if (!pcstream) return 1;
1500 if (!pTPC) return 1;
1504 static Double_t meanVertX=0;
1505 static Double_t rmsVertX=0;
1506 static Double_t meanVertY=0;
1507 static Double_t rmsVertY=0;
1508 static Double_t meanVertZ=0;
1509 static Double_t rmsVertZ=0;
1510 static Double_t vertStatus=0;
1511 static Double_t meanMult=0;
1512 static Double_t rmsMult=0;
1513 static Double_t meanMultPos=0;
1514 static Double_t rmsMultPos=0;
1515 static Double_t meanMultNeg=0;
1516 static Double_t rmsMultNeg=0;
1517 static Double_t vertAll = 0;
1518 static Double_t vertOK = 0;
1522 if (pTPC->GetHistos()->FindObject("h_tpc_event_6") && !fgForceTHnSparse) {
1523 his1D = dynamic_cast<TH1*>(pTPC->GetHistos()->FindObject("h_tpc_event_6")->Clone());
1525 his1D = pTPC->GetTPCEventHisto()->Projection(6);
1527 if(!his1D) return 1;
1529 vertAll = his1D->GetEntries();
1530 vertOK = his1D->GetBinContent(2);
1532 vertStatus = vertOK / vertAll;
1537 pTPC->GetTPCEventHisto()->GetAxis(6)->SetRange(2,2);
1539 if (pTPC->GetHistos()->FindObject("h_tpc_event_recvertex_0") && !fgForceTHnSparse) {
1540 his1D = dynamic_cast<TH1*>(pTPC->GetHistos()->FindObject("h_tpc_event_recvertex_0")->Clone());
1542 his1D = pTPC->GetTPCEventHisto()->Projection(0);
1544 if(!his1D) return 1;
1546 meanVertX = his1D->GetMean();
1547 rmsVertX = his1D->GetRMS();
1553 if (pTPC->GetHistos()->FindObject("h_tpc_event_recvertex_1") && !fgForceTHnSparse) {
1554 his1D = dynamic_cast<TH1*>(pTPC->GetHistos()->FindObject("h_tpc_event_recvertex_1")->Clone());
1556 his1D = pTPC->GetTPCEventHisto()->Projection(1);
1558 if(!his1D) return 1;
1560 meanVertY = his1D->GetMean();
1561 rmsVertY = his1D->GetRMS();
1565 if (pTPC->GetHistos()->FindObject("h_tpc_event_recvertex_2") && !fgForceTHnSparse) {
1566 hc = dynamic_cast<TH1*>(pTPC->GetHistos()->FindObject("h_tpc_event_recvertex_2"));
1568 //his1D = dynamic_cast<TH1*>(pTPC->GetHistos()->FindObject("h_tpc_event_recvertex_2")->Clone());
1569 his1D = (TH1*)hc->Clone();
1571 his1D = pTPC->GetTPCEventHisto()->Projection(2);
1573 if(!his1D) return 1;
1575 meanVertZ = his1D->GetMean();
1576 rmsVertZ = his1D->GetRMS();
1580 if (pTPC->GetHistos()->FindObject("h_tpc_event_recvertex_3") && !fgForceTHnSparse) {
1581 hc = dynamic_cast<TH1*>(pTPC->GetHistos()->FindObject("h_tpc_event_recvertex_3"));
1583 //his1D = dynamic_cast<TH1*>(pTPC->GetHistos()->FindObject("h_tpc_event_recvertex_3")->Clone());
1584 his1D = (TH1*)hc->Clone();
1586 his1D = pTPC->GetTPCEventHisto()->Projection(3);
1588 if(!his1D) return 1;
1590 meanMult = his1D->GetMean();
1591 rmsMult = his1D->GetRMS();
1595 if (pTPC->GetHistos()->FindObject("h_tpc_event_recvertex_4") && !fgForceTHnSparse) {
1596 his1D = dynamic_cast<TH1*>(pTPC->GetHistos()->FindObject("h_tpc_event_recvertex_4")->Clone());
1598 his1D = pTPC->GetTPCEventHisto()->Projection(4);
1600 if(!his1D) return 1;
1602 meanMultPos = his1D->GetMean();
1603 rmsMultPos = his1D->GetRMS();
1606 if (pTPC->GetHistos()->FindObject("h_tpc_event_recvertex_5") && !fgForceTHnSparse) {
1607 his1D = dynamic_cast<TH1*>(pTPC->GetHistos()->FindObject("h_tpc_event_recvertex_5")->Clone());
1609 his1D = pTPC->GetTPCEventHisto()->Projection(5);
1611 if(!his1D) return 1;
1613 meanMultNeg = his1D->GetMean();
1614 rmsMultNeg = his1D->GetRMS();
1617 pTPC->GetTPCEventHisto()->GetAxis(6)->SetRange(1,2);
1619 (*pcstream)<<"tpcQA"<<
1620 "meanVertX="<<meanVertX<<
1621 "rmsVertX="<<rmsVertX<<
1622 "meanVertY="<<meanVertY<<
1623 "rmsVertY="<<rmsVertY<<
1624 "meanVertZ="<<meanVertZ<<
1625 "rmsVertZ="<<rmsVertZ<<
1626 "vertStatus="<<vertStatus<<
1627 "vertAll="<<vertAll<<
1629 "meanMult="<<meanMult<<
1630 "rmsMult="<<rmsMult<<
1631 "meanMultPos="<<meanMultPos<<
1632 "rmsMultPos="<<rmsMultPos<<
1633 "meanMultNeg="<<meanMultNeg<<
1634 "rmsMultNeg="<<rmsMultNeg;
1639 //_____________________________________________________________________________
1640 Int_t AliTPCPerformanceSummary::AnalyzePt(const AliPerformanceTPC* pTPC, TTreeSRedirector* const pcstream)
1643 // Analyse DCA R imperfections for positive particles
1646 if (!pcstream) return 256;
1647 if (!pTPC) return 256;
1650 static Double_t meanPtAPos = 0;
1651 static Double_t mediumPtAPos = 0;
1652 static Double_t highPtAPos = 0;
1653 static Double_t meanPtCPos = 0;
1654 static Double_t mediumPtCPos = 0;
1655 static Double_t highPtCPos = 0;
1657 static Double_t meanPtANeg = 0;
1658 static Double_t mediumPtANeg = 0;
1659 static Double_t highPtANeg = 0;
1660 static Double_t meanPtCNeg = 0;
1661 static Double_t mediumPtCNeg = 0;
1662 static Double_t highPtCNeg = 0;
1667 if (pTPC->GetHistos()->FindObject("h_tpc_track_pos_recvertex_3_5_7")) {
1669 his3D1 = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_pos_recvertex_3_5_7"));
1670 if(!his3D1) return 256;
1672 his3D1->GetYaxis()->SetRangeUser(0.1,0.8);
1674 his3D1->GetZaxis()->SetRangeUser(0.25,10);
1675 meanPtAPos = his3D1->GetMean(3);
1676 his3D1->GetZaxis()->SetRangeUser(2,5);
1677 mediumPtAPos = his3D1->GetMean(3);
1678 his3D1->GetZaxis()->SetRangeUser(5,10);
1679 highPtAPos = his3D1->GetMean(3);
1681 his3D1->GetYaxis()->SetRangeUser(-0.8,-0.1);
1683 his3D1->GetZaxis()->SetRangeUser(0.25,10);
1684 meanPtCPos = his3D1->GetMean(3);
1685 his3D1->GetZaxis()->SetRangeUser(2,5);
1686 mediumPtCPos = his3D1->GetMean(3);
1687 his3D1->GetZaxis()->SetRangeUser(5,10);
1688 highPtCPos = his3D1->GetMean(3);
1690 his3D1->GetYaxis()->SetRangeUser(-1,1);
1691 his3D1->GetZaxis()->SetRangeUser(0.25,10);
1695 if (pTPC->GetHistos()->FindObject("h_tpc_track_neg_recvertex_3_5_7")) {
1697 his3D2 = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_neg_recvertex_3_5_7"));
1698 if(!his3D2) return 256;
1700 his3D2->GetYaxis()->SetRangeUser(0.1,0.8);
1702 his3D2->GetZaxis()->SetRangeUser(0.25,10);
1703 meanPtANeg = his3D2->GetMean(3);
1704 his3D2->GetZaxis()->SetRangeUser(2,5);
1705 mediumPtANeg = his3D2->GetMean(3);
1706 his3D2->GetZaxis()->SetRangeUser(5,10);
1707 highPtANeg = his3D2->GetMean(3);
1709 his3D2->GetYaxis()->SetRangeUser(-0.8,-0.1);
1711 his3D2->GetZaxis()->SetRangeUser(0.25,10);
1712 meanPtCNeg = his3D2->GetMean(3);
1713 his3D2->GetZaxis()->SetRangeUser(2,5);
1714 mediumPtCNeg = his3D2->GetMean(3);
1715 his3D2->GetZaxis()->SetRangeUser(5,10);
1716 highPtCNeg = his3D2->GetMean(3);
1718 his3D2->GetYaxis()->SetRangeUser(-1,1);
1719 his3D2->GetZaxis()->SetRangeUser(0.25,10);
1726 (*pcstream)<<"tpcQA"<<
1727 "meanPtAPos="<< meanPtAPos<<
1728 "mediumPtAPos="<< mediumPtAPos<<
1729 "highPtAPos="<< highPtAPos<<
1731 "meanPtCPos="<< meanPtCPos<<
1732 "mediumPtCPos="<< mediumPtCPos<<
1733 "highPtCPos="<< highPtCPos<<
1735 "meanPtANeg="<< meanPtANeg<<
1736 "mediumPtANeg="<< mediumPtANeg<<
1737 "highPtANeg="<< highPtANeg<<
1739 "meanPtCNeg="<< meanPtCNeg<<
1740 "mediumPtCNeg="<< mediumPtCNeg<<
1741 "highPtCNeg="<< highPtCNeg;
1747 //_____________________________________________________________________________
1749 Int_t AliTPCPerformanceSummary::AnalyzeChargeOverPt(const AliPerformanceTPC* pTPC, TTreeSRedirector* const pcstream){
1751 // Analyse DCA R imperfections for positive particles
1754 if (!pcstream) return 512;
1755 if (!pTPC) return 512;
1758 static Double_t qOverPt = 0;
1759 static Double_t qOverPtA = 0;
1760 static Double_t qOverPtC = 0;
1766 TF1 *fp1 = new TF1("fp1","pol2",-1.0,1.0);
1767 TF1 *fp2 = new TF1("fp2","pol2",-1.0,1.0);
1768 TF1 *fp3 = new TF1("fp3","pol2",-1.0,1.0);
1770 if (pTPC->GetHistos()->FindObject("h_tpc_track_all_recvertex_5_8")) {
1772 his2D = dynamic_cast<TH2*>(pTPC->GetHistos()->FindObject("h_tpc_track_all_recvertex_5_8"));
1773 if(!his2D) return 512;
1775 his1D1 = his2D->ProjectionX();
1776 his1D1->Fit(fp1,"R");
1777 if(fp1->GetParameter(2)!=0){
1778 qOverPt = (-1.0)*(fp1->GetParameter(1)/(2.0*fp1->GetParameter(2)));
1783 his2D->GetYaxis()->SetRangeUser(0.1,0.8);
1784 his1D2 = his2D->ProjectionX();
1785 his1D2->Fit(fp2,"R");
1786 if(fp2->GetParameter(2)!=0)
1787 qOverPtA = (-1.0)*(fp2->GetParameter(1)/(2.0*fp2->GetParameter(2)));
1791 his2D->GetYaxis()->SetRangeUser(-0.8,-0.1);
1792 his1D3 = his2D->ProjectionX();
1793 his1D3->Fit(fp3,"R");
1794 if(fp3->GetParameter(2)!=0)
1795 qOverPtC = (-1.0)*(fp3->GetParameter(1)/(2.0*fp3->GetParameter(2)));
1799 his2D->GetYaxis()->SetRangeUser(-1.0,1.0);
1803 (*pcstream)<<"tpcQA"<<
1804 "qOverPt="<< qOverPt<<
1805 "qOverPtA="<< qOverPtA<<
1806 "qOverPtC="<< qOverPtC;
1811 Int_t AliTPCPerformanceSummary::AnalyzeMatch(const AliPerformanceMatch* pMatch, TTreeSRedirector* const pcstream)
1813 /* if ((pMatch == 0) or (0 == pcstream)) { printf("this will not work anyway..."); }
1814 printf("funtion not implemented");*/
1816 if (!pcstream) return 1024;
1817 if (!pMatch) return 1024;
1818 static Double_t tpcItsMatchA = 0;
1819 static Double_t tpcItsMatchHighPtA = 0;
1820 static Double_t tpcItsMatchC = 0;
1821 static Double_t tpcItsMatchHighPtC = 0;
1825 if(pMatch->GetHistos()->FindObject("h_tpc_match_trackingeff_all_2_3") &&
1826 pMatch->GetHistos()->FindObject("h_tpc_match_trackingeff_tpc_2_3")){
1827 h2D = dynamic_cast<TH2*>(pMatch->GetHistos()->FindObject("h_tpc_match_trackingeff_all_2_3"));
1828 h2D1 = dynamic_cast<TH2*>(pMatch->GetHistos()->FindObject("h_tpc_match_trackingeff_tpc_2_3"));
1833 h2D->GetXaxis()->SetRangeUser(0,1.5);
1834 h2D1->GetXaxis()->SetRangeUser(0,1.5);
1836 Double_t entries,entries1;
1837 entries = h2D->GetEffectiveEntries();
1838 entries1 = h2D1->GetEffectiveEntries();
1840 tpcItsMatchA = entries1/entries;
1842 h2D->GetYaxis()->SetRangeUser(4.01,20.);
1843 h2D1->GetYaxis()->SetRangeUser(4.01,20.);
1844 entries = h2D->GetEffectiveEntries();
1845 entries1 = h2D1->GetEffectiveEntries();
1847 tpcItsMatchHighPtA = entries1/entries;
1850 h2D->GetXaxis()->SetRangeUser(-1.5,-0.01);
1851 h2D1->GetXaxis()->SetRangeUser(-1.5,-0.01);
1852 h2D->GetYaxis()->SetRangeUser(0.0,20.);
1853 h2D1->GetYaxis()->SetRangeUser(0.0,20.);
1855 entries = h2D->GetEffectiveEntries();
1856 entries1 = h2D1->GetEffectiveEntries();
1858 tpcItsMatchC = entries1/entries;
1860 h2D->GetYaxis()->SetRangeUser(4.01,20.);
1861 h2D1->GetYaxis()->SetRangeUser(4.01,20.);
1862 entries = h2D->GetEffectiveEntries();
1863 entries1 = h2D1->GetEffectiveEntries();
1865 tpcItsMatchHighPtC = entries1/entries;
1867 h2D->GetXaxis()->SetRangeUser(-1.5,1.5);
1868 h2D1->GetXaxis()->SetRangeUser(-1.5,1.5);
1869 h2D->GetYaxis()->SetRangeUser(0.0,20.);
1870 h2D1->GetYaxis()->SetRangeUser(0.0,20.);
1875 (*pcstream)<<"tpcQA"<<
1876 "tpcItsMatchA="<< tpcItsMatchA<<
1877 "tpcItsMatchHighPtA="<< tpcItsMatchHighPtA<<
1878 "tpcItsMatchC="<< tpcItsMatchC<<
1879 "tpcItsMatchHighPtC="<< tpcItsMatchHighPtC;
1884 Int_t AliTPCPerformanceSummary::AnalyzePull(const AliPerformanceMatch* pPull, TTreeSRedirector* const pcstream)
1886 /* if ((pPull == 0) or (0 == pcstream)) { printf("this will not work anyway..."); }
1887 printf("funtion not implemented");*/
1889 if (!pcstream) return 2048;
1890 if (!pPull) return 2048;
1891 static Double_t phiPull = 0;
1892 static Double_t phiPullHighPt = 0;
1893 static Double_t ptPull = 0;
1894 static Double_t ptPullHighPt = 0;
1895 static Double_t yPull = 0;
1896 static Double_t yPullHighPt = 0;
1897 static Double_t zPull = 0;
1898 static Double_t zPullHighPt = 0;
1899 static Double_t lambdaPull = 0;
1900 static Double_t lambdaPullHighPt = 0;
1903 if(pPull->GetHistos()->FindObject("h_tpc_match_pull_2_7")){
1904 h2D1 = dynamic_cast<TH2*>(pPull->GetHistos()->FindObject("h_tpc_match_pull_2_7"));
1906 phiPull = h2D1->GetMean(2);
1907 h2D1->SetAxisRange(0.0,1.0/5.0,"X");
1908 phiPullHighPt = h2D1->GetMean(2);
1909 h2D1->SetAxisRange(0.0,10.0,"X");
1914 if(pPull->GetHistos()->FindObject("h_tpc_match_pull_4_7")){
1915 h2D2 = dynamic_cast<TH2*>(pPull->GetHistos()->FindObject("h_tpc_match_pull_4_7"));
1917 ptPull = h2D2->GetMean(2);
1919 h2D2->SetAxisRange(0.0,1.0/5.0,"X");
1920 ptPullHighPt = h2D2->GetMean(2);
1921 h2D2->SetAxisRange(0.0,10.0,"X");
1926 if(pPull->GetHistos()->FindObject("h_tpc_match_pull_0_7")){
1927 h2D3 = dynamic_cast<TH2*>(pPull->GetHistos()->FindObject("h_tpc_match_pull_0_7"));
1929 yPull = h2D3->GetMean(2);
1931 h2D3->SetAxisRange(0.0,1.0/5.0,"X");
1932 yPullHighPt = h2D3->GetMean(2);
1933 h2D3->SetAxisRange(0.0,10.0,"X");
1938 if(pPull->GetHistos()->FindObject("h_tpc_match_pull_1_7")){
1939 h2D4 = dynamic_cast<TH2*>(pPull->GetHistos()->FindObject("h_tpc_match_pull_1_7"));
1941 zPull = h2D4->GetMean(2);
1943 h2D4->SetAxisRange(0.0,1.0/5.0,"X");
1944 zPullHighPt = h2D4->GetMean(2);
1945 h2D4->SetAxisRange(0.0,10.0,"X");
1950 if(pPull->GetHistos()->FindObject("h_tpc_match_pull_3_7")){
1951 h2D5 = dynamic_cast<TH2*>(pPull->GetHistos()->FindObject("h_tpc_match_pull_3_7"));
1953 lambdaPull = h2D5->GetMean(2);
1955 h2D5->SetAxisRange(0.0,1.0/5.0,"X");
1956 lambdaPullHighPt = h2D5->GetMean(2);
1957 h2D5->SetAxisRange(0.0,10.0,"X");
1961 (*pcstream)<<"tpcQA"<<
1962 "phiPull="<< phiPull<<
1963 "phiPullHighPt="<< phiPullHighPt<<
1964 "ptPull="<< ptPull<<
1965 "ptPullHighPt="<< ptPullHighPt<<
1967 "yPullHighPt="<< yPullHighPt<<
1969 "zPullHighPt="<< zPullHighPt<<
1970 "lambdaPull="<< lambdaPull<<
1971 "lambdaPullHighPt="<< lambdaPullHighPt;
1975 Int_t AliTPCPerformanceSummary::AnalyzeConstrain(const AliPerformanceMatch* pConstrain, TTreeSRedirector* pcstream)
1977 if (!pcstream) return 5126;
1978 if (!pConstrain) return 5126;
1981 static Double_t tpcConstrainPhiA = 0;
1982 static Double_t tpcConstrainPhiC = 0;
1984 if (pConstrain->GetHistos()->FindObject("h_tpc_constrain_tpc_0_2_3")) {
1986 his3D = dynamic_cast<TH3*>(pConstrain->GetHistos()->FindObject("h_tpc_constrain_tpc_0_2_3"));//phi pull:pt:eta
1987 if(!his3D) return 5126;
1989 his3D->GetZaxis()->SetRangeUser(0.0,1.0);
1990 tpcConstrainPhiA = his3D->GetMean(1);
1991 his3D->GetZaxis()->SetRangeUser(-1.0,-0.001);
1992 tpcConstrainPhiC = his3D->GetMean(1);
1995 (*pcstream)<<"tpcQA"<<
1996 "tpcConstrainPhiA="<<tpcConstrainPhiA <<
1997 "tpcConstrainPhiC="<< tpcConstrainPhiC;
2002 //_____________________________________________________________________________
2003 Int_t AliTPCPerformanceSummary::AnalyzeQAPosNegDpT(const AliPerformanceTPC* pTPC, TTreeSRedirector* const pcstream)
2005 //function which plot 1/Pt for negative and
2006 //positive particles
2008 if (!pcstream) return 512;
2009 if (!pTPC) return 512;
2019 static Double_t deltaPtC = 0;
2020 static Double_t deltaPtchi2C = 0;
2021 static Double_t slopeC = 0;
2022 static Double_t deltaPtA = 0;
2023 static Double_t deltaPtchi2A = 0;
2024 static Double_t slopeA = 0;
2025 static Double_t deltaPt = 0;
2026 static Double_t deltaPtchi2 = 0;
2027 static Double_t slope = 0;
2028 static Double_t deltaPt_Err = 0;
2029 static Double_t deltaPtA_Err = 0;
2030 static Double_t deltaPtC_Err = 0;
2035 if(pTPC->GetHistos()->FindObject("h_tpc_track_pos_recvertex_0_5_7"))
2037 pos3 = dynamic_cast<TH3D*>(pTPC->GetHistos()->FindObject("h_tpc_track_pos_recvertex_0_5_7"));
2038 if(!pos3) return 512;
2040 pos = pos3->ProjectionZ("pos",71,-1,6,25);
2041 posC = pos3->ProjectionZ("posC",71,-1,6,15);
2042 posA = pos3->ProjectionZ("posA",71,-1,16,25);
2045 if(pTPC->GetHistos()->FindObject("h_tpc_track_neg_recvertex_0_5_7")){
2046 neg3 = dynamic_cast<TH3D*>(pTPC->GetHistos()->FindObject("h_tpc_track_neg_recvertex_0_5_7"));
2047 if(!neg3) return 512;
2049 neg = neg3->ProjectionZ("neg",71,-1,6,25);
2050 negC = neg3->ProjectionZ("negC",71,-1,6,15);
2051 negA = neg3->ProjectionZ("negA",71,-1,16,25);
2062 pos->Scale(1.,"width");
2063 neg->Scale(1.,"width");
2064 posA->Scale(1.,"width");
2065 negA->Scale(1.,"width");
2066 posC->Scale(1.,"width");
2067 negC->Scale(1.,"width");
2071 TF1 fpt("fpt","[1]*exp(-1/((1/x))*[0])",0.1,10);
2072 TF1 fpt2("fpt2","[1]*exp(-1/((1/x))*[0])",0.1,10);
2073 fpt.SetParameters(1,0.5);
2074 fpt2.SetParameters(1,0.5);
2075 pos->Fit(&fpt,"","",1,4); pos->Fit(&fpt,"","",1,4); pos->Fit(&fpt,"","",1,4);
2076 neg->Fit(&fpt2,"","",1,4); neg->Fit(&fpt2,"","",1,4); neg->Fit(&fpt2,"","",1,4);
2078 slope = (fpt.GetParameter(0)+fpt2.GetParameter(0))/2.;
2080 TH1D* ratio = new TH1D(*pos);
2084 TF1 fptRatio("fptratio","[2]*exp(-1/((1/x)+[1])*[0])/exp(-1/((1/x)-[1])*[0])",0.1,10);
2085 fptRatio.SetParameters(0.5,0.006,1);
2086 fptRatio.FixParameter(0,slope);
2088 ratio->Fit(&fptRatio,"","",1,4); ratio->Fit(&fptRatio,"","",1,4);
2089 ratio->Fit(&fptRatio,"","",1,4);
2091 deltaPt = fptRatio.GetParameter(1);
2092 deltaPtchi2 = fptRatio.GetChisquare();
2095 deltaPt_Err = fptRatio.GetParError(1);
2100 TF1 fptA("fptA","[1]*exp(-1/((1/x))*[0])",0.1,10);
2101 TF1 fpt2A("fpt2A","[1]*exp(-1/((1/x))*[0])",0.1,10);
2102 fptA.SetParameters(1,0.5);
2103 fpt2A.SetParameters(1,0.5);
2104 posA->Fit(&fptA,"","",1,4); posA->Fit(&fptA,"","",1,4); posA->Fit(&fptA,"","",1,4);
2105 negA->Fit(&fpt2A,"","",1,4); negA->Fit(&fpt2A,"","",1,4); negA->Fit(&fpt2A,"","",1,4);
2107 slopeA = (fptA.GetParameter(0)+fpt2A.GetParameter(0))/2.;
2109 TH1D* ratioA = new TH1D(*posA);
2110 ratioA->Divide(negA);
2113 TF1 fptRatioA("fptratioA","[2]*exp(-1/((1/x)+[1])*[0])/exp(-1/((1/x)-[1])*[0])",0.1,10);
2114 fptRatioA.SetParameters(0.5,0.006,1);
2115 fptRatioA.FixParameter(0,slopeA);
2117 ratioA->Fit(&fptRatioA,"","",1,4); ratio->Fit(&fptRatioA,"","",1,4);
2118 ratioA->Fit(&fptRatioA,"","",1,4);
2120 deltaPtA = fptRatioA.GetParameter(1);
2121 deltaPtchi2A = fptRatioA.GetChisquare();
2124 deltaPtA_Err = fptRatioA.GetParError(1);
2132 TF1 fptC("fptC","[1]*exp(-1/((1/x))*[0])",0.1,10);
2133 TF1 fpt2C("fpt2C","[1]*exp(-1/((1/x))*[0])",0.1,10);
2134 fptC.SetParameters(1,0.5);
2135 fpt2C.SetParameters(1,0.5);
2136 posC->Fit(&fptC,"","",1,4); posC->Fit(&fptC,"","",1,4); posC->Fit(&fptC,"","",1,4);
2137 negC->Fit(&fpt2C,"","",1,4); negC->Fit(&fpt2C,"","",1,4); negC->Fit(&fpt2C,"","",1,4);
2139 slopeC = (fptC.GetParameter(0)+fpt2C.GetParameter(0))/2.;
2141 TH1D* ratioC = new TH1D(*posC);
2142 ratioC->Divide(negC);
2145 TF1 fptRatioC("fptratioC","[2]*exp(-1/((1/x)+[1])*[0])/exp(-1/((1/x)-[1])*[0])",0.1,10);
2146 fptRatioC.SetParameters(0.5,0.006,1);
2147 fptRatioC.FixParameter(0,slopeC);
2149 ratioC->Fit(&fptRatioC,"","",1,4); ratio->Fit(&fptRatioC,"","",1,4);
2150 ratioC->Fit(&fptRatioC,"","",1,4);
2152 deltaPtC = fptRatioC.GetParameter(1);
2153 deltaPtchi2C = fptRatioC.GetChisquare();
2156 deltaPtC_Err = fptRatioC.GetParError(1);
2163 (*pcstream)<<"tpcQA"<<
2164 "deltaPt="<< deltaPt<<
2165 "deltaPtchi2="<< deltaPtchi2<<
2166 "deltaPtA="<< deltaPtA<<
2167 "deltaPtchi2A="<< deltaPtchi2A<<
2168 "deltaPtC="<< deltaPtC<<
2169 "deltaPtchi2C="<< deltaPtchi2C<<
2170 "deltaPt_Err="<< deltaPt_Err<<
2171 "deltaPtA_Err="<< deltaPtA_Err<<
2172 "deltaPtC_Err="<< deltaPtC_Err;
2177 //_____________________________________________________________________________
2178 Int_t AliTPCPerformanceSummary::AnalyzeQADCAFitParameter(const AliPerformanceTPC* pTPC, TTreeSRedirector* const pcstream)
2182 //function which retrieve DCA fit parameters
2185 if (!pcstream) return 16;
2186 if (!pTPC) return 16;
2193 static Double_t dcar_posA_0=0;
2194 static Double_t dcar_posA_1=0;
2195 static Double_t dcar_posA_2=0;
2196 static Double_t dcar_posA_chi2=0;
2197 static Double_t dcar_posA_0_Err=0;
2198 static Double_t dcar_posA_1_Err=0;
2199 static Double_t dcar_posA_2_Err=0;
2201 static Double_t dcar_posC_0=0;
2202 static Double_t dcar_posC_1=0;
2203 static Double_t dcar_posC_2=0;
2204 static Double_t dcar_posC_chi2=0;
2205 static Double_t dcar_posC_0_Err=0;
2206 static Double_t dcar_posC_1_Err=0;
2207 static Double_t dcar_posC_2_Err=0;
2209 static Double_t dcaz_posA_0=0;
2210 static Double_t dcaz_posA_1=0;
2211 static Double_t dcaz_posA_2=0;
2212 static Double_t dcaz_posA_chi2=0;
2213 static Double_t dcaz_posA_0_Err=0;
2214 static Double_t dcaz_posA_1_Err=0;
2215 static Double_t dcaz_posA_2_Err=0;
2217 static Double_t dcaz_posC_0=0;
2218 static Double_t dcaz_posC_1=0;
2219 static Double_t dcaz_posC_2=0;
2220 static Double_t dcaz_posC_chi2=0;
2221 static Double_t dcaz_posC_0_Err=0;
2222 static Double_t dcaz_posC_1_Err=0;
2223 static Double_t dcaz_posC_2_Err=0;
2225 static Double_t dcar_negA_0=0;
2226 static Double_t dcar_negA_1=0;
2227 static Double_t dcar_negA_2=0;
2228 static Double_t dcar_negA_chi2=0;
2229 static Double_t dcar_negA_0_Err=0;
2230 static Double_t dcar_negA_1_Err=0;
2231 static Double_t dcar_negA_2_Err=0;
2233 static Double_t dcar_negC_0=0;
2234 static Double_t dcar_negC_1=0;
2235 static Double_t dcar_negC_2=0;
2236 static Double_t dcar_negC_chi2=0;
2237 static Double_t dcar_negC_0_Err=0;
2238 static Double_t dcar_negC_1_Err=0;
2239 static Double_t dcar_negC_2_Err=0;
2241 static Double_t dcaz_negA_0=0;
2242 static Double_t dcaz_negA_1=0;
2243 static Double_t dcaz_negA_2=0;
2244 static Double_t dcaz_negA_chi2=0;
2245 static Double_t dcaz_negA_0_Err=0;
2246 static Double_t dcaz_negA_1_Err=0;
2247 static Double_t dcaz_negA_2_Err=0;
2249 static Double_t dcaz_negC_0=0;
2250 static Double_t dcaz_negC_1=0;
2251 static Double_t dcaz_negC_2=0;
2252 static Double_t dcaz_negC_chi2=0;
2253 static Double_t dcaz_negC_0_Err=0;
2254 static Double_t dcaz_negC_1_Err=0;
2255 static Double_t dcaz_negC_2_Err=0;
2257 if (pTPC->GetHistos()->FindObject("h_tpc_track_pos_recvertex_3_5_6")) {
2258 dcar_pos3 = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_pos_recvertex_3_5_6"));
2261 if (pTPC->GetHistos()->FindObject("h_tpc_track_pos_recvertex_3_5_6")) {
2262 dcaz_pos3 = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_pos_recvertex_4_5_6"));
2265 if (pTPC->GetHistos()->FindObject("h_tpc_track_pos_recvertex_3_5_6")) {
2266 dcar_neg3 = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_neg_recvertex_3_5_6"));
2269 if (pTPC->GetHistos()->FindObject("h_tpc_track_pos_recvertex_3_5_6")) {
2270 dcaz_neg3 = dynamic_cast<TH3*>(pTPC->GetHistos()->FindObject("h_tpc_track_neg_recvertex_4_5_6"));
2273 TF1 fit("fit","[0]+[1]*cos(x)+[2]*sin(x)",0,7);
2275 dcar_pos3->GetYaxis()->SetRangeUser(0,0.99);
2276 TH1* dcar_posA = (dynamic_cast<TH2*>(dcar_pos3->Project3D("xz_1")))->ProfileX();
2277 dcar_posA->Fit(&fit,"NQ");
2278 dcar_posA_0 = fit.GetParameter(0);
2279 dcar_posA_1 = fit.GetParameter(1);
2280 dcar_posA_2 = fit.GetParameter(2);
2281 dcar_posA_chi2 = fit.GetChisquare();
2282 dcar_posA_0_Err = fit.GetParError(0);
2283 dcar_posA_1_Err = fit.GetParError(1);
2284 dcar_posA_2_Err = fit.GetParError(2);
2286 dcar_pos3->GetYaxis()->SetRangeUser(-1.0,-0.01);
2287 TH1* dcar_posC = (dynamic_cast<TH2*>(dcar_pos3->Project3D("xz_2")))->ProfileX();
2288 dcar_posC->Fit(&fit,"NQ");
2289 dcar_posC_0 = fit.GetParameter(0);
2290 dcar_posC_1 = fit.GetParameter(1);
2291 dcar_posC_2 = fit.GetParameter(2);
2292 dcar_posC_chi2 = fit.GetChisquare();
2293 dcar_posC_0_Err = fit.GetParError(0);
2294 dcar_posC_1_Err = fit.GetParError(1);
2295 dcar_posC_2_Err = fit.GetParError(2);
2297 dcaz_pos3->GetYaxis()->SetRangeUser(0,0.99);
2298 TH1* dcaz_posA = (dynamic_cast<TH2*>(dcaz_pos3->Project3D("xz_3")))->ProfileX();
2299 dcaz_posA->Fit(&fit,"NQ");
2300 dcaz_posA_0 = fit.GetParameter(0);
2301 dcaz_posA_1 = fit.GetParameter(1);
2302 dcaz_posA_2 = fit.GetParameter(2);
2303 dcaz_posA_chi2 = fit.GetChisquare();
2304 dcaz_posA_0_Err = fit.GetParError(0);
2305 dcaz_posA_1_Err = fit.GetParError(1);
2306 dcaz_posA_2_Err = fit.GetParError(2);
2308 dcaz_pos3->GetYaxis()->SetRangeUser(-1.0,-0.01);
2309 TH1* dcaz_posC = (dynamic_cast<TH2*>(dcaz_pos3->Project3D("xz_4")))->ProfileX();
2310 dcaz_posC->Fit(&fit,"NQ");
2311 dcaz_posC_0 = fit.GetParameter(0);
2312 dcaz_posC_1 = fit.GetParameter(1);
2313 dcaz_posC_2 = fit.GetParameter(2);
2314 dcaz_posC_chi2 = fit.GetChisquare();
2315 dcaz_posC_0_Err = fit.GetParError(0);
2316 dcaz_posC_1_Err = fit.GetParError(1);
2317 dcaz_posC_2_Err = fit.GetParError(2);
2321 dcar_neg3->GetYaxis()->SetRangeUser(0,0.99);
2322 TH1* dcar_negA = (dynamic_cast<TH2*>(dcar_neg3->Project3D("xz_1")))->ProfileX();
2323 dcar_negA->Fit(&fit,"NQ");
2324 dcar_negA_0 = fit.GetParameter(0);
2325 dcar_negA_1 = fit.GetParameter(1);
2326 dcar_negA_2 = fit.GetParameter(2);
2327 dcar_negA_chi2 = fit.GetChisquare();
2328 dcar_negA_0_Err = fit.GetParError(0);
2329 dcar_negA_1_Err = fit.GetParError(1);
2330 dcar_negA_2_Err = fit.GetParError(2);
2332 dcar_neg3->GetYaxis()->SetRangeUser(-1.0,-0.01);
2333 TH1* dcar_negC = (dynamic_cast<TH2*>(dcar_neg3->Project3D("xz_2")))->ProfileX();
2334 dcar_negC->Fit(&fit,"NQ");
2335 dcar_negC_0 = fit.GetParameter(0);
2336 dcar_negC_1 = fit.GetParameter(1);
2337 dcar_negC_2 = fit.GetParameter(2);
2338 dcar_negC_chi2 = fit.GetChisquare();
2339 dcar_negC_0_Err = fit.GetParError(0);
2340 dcar_negC_1_Err = fit.GetParError(1);
2341 dcar_negC_2_Err = fit.GetParError(2);
2343 dcaz_neg3->GetYaxis()->SetRangeUser(0,0.99);
2344 TH1* dcaz_negA = (dynamic_cast<TH2*>(dcaz_neg3->Project3D("xz_3")))->ProfileX();
2345 dcaz_negA->Fit(&fit,"NQ");
2346 dcaz_negA_0 = fit.GetParameter(0);
2347 dcaz_negA_1 = fit.GetParameter(1);
2348 dcaz_negA_2 = fit.GetParameter(2);
2349 dcaz_negA_chi2 = fit.GetChisquare();
2350 dcaz_negA_0_Err = fit.GetParError(0);
2351 dcaz_negA_1_Err = fit.GetParError(1);
2352 dcaz_negA_2_Err = fit.GetParError(2);
2354 dcaz_neg3->GetYaxis()->SetRangeUser(-1.0,-0.01);
2355 TH1* dcaz_negC = (dynamic_cast<TH2*>(dcaz_neg3->Project3D("xz_4")))->ProfileX();
2356 dcaz_negC->Fit(&fit,"NQ");
2357 dcaz_negC_0 = fit.GetParameter(0);
2358 dcaz_negC_1 = fit.GetParameter(1);
2359 dcaz_negC_2 = fit.GetParameter(2);
2360 dcaz_negC_chi2 = fit.GetChisquare();
2361 dcaz_negC_0_Err = fit.GetParError(0);
2362 dcaz_negC_1_Err = fit.GetParError(1);
2363 dcaz_negC_2_Err = fit.GetParError(2);
2366 // store results (shift in dca) in ttree
2368 (*pcstream)<<"tpcQA"<<
2369 "dcar_posA_0="<< dcar_posA_0<<
2370 "dcar_posA_1="<< dcar_posA_1<<
2371 "dcar_posA_2="<< dcar_posA_2<<
2372 "dcar_posA_chi2="<< dcar_posA_chi2<<
2373 "dcar_posA_0_Err="<< dcar_posA_0_Err<<
2374 "dcar_posA_1_Err="<< dcar_posA_1_Err<<
2375 "dcar_posA_2_Err="<< dcar_posA_2_Err;
2377 (*pcstream)<<"tpcQA"<<
2378 "dcaz_posA_0="<< dcaz_posA_0<<
2379 "dcaz_posA_1="<< dcaz_posA_1<<
2380 "dcaz_posA_2="<< dcaz_posA_2<<
2381 "dcaz_posA_chi2="<< dcaz_posA_chi2<<
2382 "dcaz_posA_0_Err="<< dcaz_posA_0_Err<<
2383 "dcaz_posA_1_Err="<< dcaz_posA_1_Err<<
2384 "dcaz_posA_2_Err="<< dcaz_posA_2_Err;
2386 (*pcstream)<<"tpcQA"<<
2387 "dcaz_posC_0="<< dcaz_posC_0<<
2388 "dcaz_posC_1="<< dcaz_posC_1<<
2389 "dcaz_posC_2="<< dcaz_posC_2<<
2390 "dcaz_posC_chi2="<< dcaz_posC_chi2<<
2391 "dcaz_posC_0_Err="<< dcaz_posC_0_Err<<
2392 "dcaz_posC_1_Err="<< dcaz_posC_1_Err<<
2393 "dcaz_posC_2_Err="<< dcaz_posC_2_Err;
2395 (*pcstream)<<"tpcQA"<<
2396 "dcar_posC_0="<< dcar_posC_0<<
2397 "dcar_posC_1="<< dcar_posC_1<<
2398 "dcar_posC_2="<< dcar_posC_2<<
2399 "dcar_posC_chi2="<< dcar_posC_chi2<<
2400 "dcar_posC_0_Err="<< dcar_posC_0_Err<<
2401 "dcar_posC_1_Err="<< dcar_posC_1_Err<<
2402 "dcar_posC_2_Err="<< dcar_posC_2_Err;
2405 (*pcstream)<<"tpcQA"<<
2406 "dcar_negA_0="<< dcar_negA_0<<
2407 "dcar_negA_1="<< dcar_negA_1<<
2408 "dcar_negA_2="<< dcar_negA_2<<
2409 "dcar_negA_chi2="<< dcar_negA_chi2<<
2410 "dcar_negA_0_Err="<< dcar_negA_0_Err<<
2411 "dcar_negA_1_Err="<< dcar_negA_1_Err<<
2412 "dcar_negA_2_Err="<< dcar_negA_2_Err;
2414 (*pcstream)<<"tpcQA"<<
2415 "dcaz_negA_0="<< dcaz_negA_0<<
2416 "dcaz_negA_1="<< dcaz_negA_1<<
2417 "dcaz_negA_2="<< dcaz_negA_2<<
2418 "dcaz_negA_chi2="<< dcaz_negA_chi2<<
2419 "dcaz_negA_0_Err="<< dcaz_negA_0_Err<<
2420 "dcaz_negA_1_Err="<< dcaz_negA_1_Err<<
2421 "dcaz_negA_2_Err="<< dcaz_negA_2_Err;
2423 (*pcstream)<<"tpcQA"<<
2424 "dcaz_negC_0="<< dcaz_negC_0<<
2425 "dcaz_negC_1="<< dcaz_negC_1<<
2426 "dcaz_negC_2="<< dcaz_negC_2<<
2427 "dcaz_negC_chi2="<< dcaz_negC_chi2<<
2428 "dcaz_negC_0_Err="<< dcaz_negC_0_Err<<
2429 "dcaz_negC_1_Err="<< dcaz_negC_1_Err<<
2430 "dcaz_negC_2_Err="<< dcaz_negC_2_Err;
2432 (*pcstream)<<"tpcQA"<<
2433 "dcar_negC_0="<< dcar_negC_0<<
2434 "dcar_negC_1="<< dcar_negC_1<<
2435 "dcar_negC_2="<< dcar_negC_2<<
2436 "dcar_negC_chi2="<< dcar_negC_chi2<<
2437 "dcar_negC_0_Err="<< dcar_negC_0_Err<<
2438 "dcar_negC_1_Err="<< dcar_negC_1_Err<<
2439 "dcar_negC_2_Err="<< dcar_negC_2_Err;