]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCcalibSummary.cxx
AliTPCcalibDButil and Summary text:
[u/mrichter/AliRoot.git] / TPC / AliTPCcalibSummary.cxx
1 /*
2 // Make a summary information of calibration.
3 // Store results in the summary trees
4 // OCDB configuration
5
6 Example usage:
7
8 gSystem->Load("libANALYSIS");
9 gSystem->Load("libTPCcalib");
10
11 Int_t irun=119037;
12 gROOT->LoadMacro("$ALICE_ROOT/TPC/scripts/OCDBscan/ConfigOCDB.C");
13 ConfigOCDB(irun)
14
15 AliTPCcalibSummary *calibSummary = new AliTPCcalibSummary;
16 calibSummary->ProcessRun(irun);
17 delete calibSummary;
18
19 */
20
21 #include <TROOT.h>
22 #include <iostream>
23 #include <fstream>
24 #include <stdio.h>
25 #include <AliCDBManager.h>
26 #include <AliCDBEntry.h>
27 #include <AliLog.h>
28 #include <AliMagF.h>
29 #include "AliTPCcalibDB.h"
30 #include "AliTPCcalibDButil.h"
31 #include "AliTPCAltroMapping.h"
32 #include "AliTPCExB.h"
33 #include "AliTPCCalROC.h"
34 #include "AliTPCCalPad.h"
35 #include "AliTPCSensorTempArray.h"
36 #include "AliGRPObject.h"
37 #include "AliTPCTransform.h"
38 #include "TFile.h"
39 #include "TKey.h"
40 #include "TObjArray.h"
41 #include "TObjString.h"
42 #include "TString.h"
43 #include "AliTPCCalPad.h"
44 #include "AliTPCROC.h"
45 #include "AliTPCParam.h"
46 #include "AliTPCCalibPulser.h"
47 #include "AliTPCCalibPedestal.h"
48 #include "AliTPCCalibCE.h"
49 #include "AliTPCExBFirst.h"
50 #include "TTreeStream.h"
51 #include "AliTPCTempMap.h"
52 #include "TVectorD.h"
53 #include "TMatrixD.h"
54 #include "AliTPCCalibRaw.h"
55 #include "AliSplineFit.h"
56 #include "TGraphErrors.h"
57 #include <AliCTPTimeParams.h>
58 #include <AliTPCcalibSummary.h>
59 #include <TStatToolkit.h>
60 #include <TCut.h>
61 //
62 //
63 //
64 AliTPCcalibSummary::AliTPCcalibSummary():
65   TNamed(),
66   fCalibDB(0),
67   fDButil(0),
68   fPcstream(0)
69 {
70   //
71   // default constructor
72   // OCDB have to be setupe before - not part of class
73   // usualy ConfigOCDB.C macro used
74   // 
75   fPcstream = new TTreeSRedirector("dcsTime.root");
76   fCalibDB = AliTPCcalibDB::Instance();
77   fDButil= new AliTPCcalibDButil;
78 }
79
80 AliTPCcalibSummary::~AliTPCcalibSummary(){
81   //
82   // destructor  - close streamer
83   //
84   delete fPcstream;  
85 }
86
87 void AliTPCcalibSummary::Process(const char * runList, Int_t first, Int_t last){
88   //
89   // runList - listOfRuns to process
90   // first   - first run to process
91   // last    - last  to process
92   // 
93   //
94   // make list of runs
95   //
96
97   ifstream inputFile;
98   inputFile.open("run.list");
99   Int_t irun=0;
100   TArrayI runArray(100000);
101   Int_t indexes[100000];
102   Int_t nruns=0;
103   printf("Runs to process:\n");
104   if (!inputFile.is_open()) { 
105     printf("Problem to open file %s\n",runList);
106   }
107   while( inputFile.good() ) {
108     inputFile >> irun;
109     printf("Run \t%d\n",irun);
110     if (irun<first) continue;  // process only subset of list
111     if (last>0 && irun>last) continue;  // process only subset of list
112     runArray[nruns]=irun;
113     nruns++;
114   }
115
116
117   TMath::Sort(nruns, runArray.fArray, indexes,kFALSE);
118   Int_t startTime = 0;
119   Int_t endTime   = 0;
120   for (Int_t run=0; run<nruns; run++){
121     irun=runArray[indexes[run]];
122     printf("Processing run %d ...\n",irun);
123     fCalibDB->SetRun(irun);
124     fDButil->UpdateFromCalibDB();
125     fDButil->SetReferenceRun(irun);
126     fDButil->UpdateRefDataFromOCDB();
127     fCalibDB->CreateGUITree("calPads.root");
128     fDButil->CreateGUIRefTree("calPadsRef.root");
129     //
130     AliDCSSensorArray *arrHV=fCalibDB->GetVoltageSensors(irun);
131     if (!arrHV) continue;
132     for  (Int_t isenHV=0; isenHV<arrHV->NumSensors(); ++isenHV){
133       AliDCSSensor *senHV=arrHV->GetSensorNum(isenHV);
134       if (!senHV) {
135         printf("Not interesting OCDB info\n");
136         continue;
137       }
138       startTime=senHV->GetStartTime();
139       endTime  =senHV->GetEndTime();
140       if (startTime>0&&endTime>0) break;
141     }
142     AliDCSSensorArray* goofieArray = fCalibDB->GetGoofieSensors(irun);  
143     if (goofieArray) fDButil->FilterGoofie(goofieArray,0.5,4.,6.85,7.05,fPcstream);
144     // don't filter goofie for the moment
145     ProcessRun(irun, startTime,endTime);
146   }
147 }
148
149
150 void AliTPCcalibSummary::ProcessRun(Int_t irun, Int_t startTime, Int_t endTime){
151   //
152   // Process run irun
153   // 
154   fCalibDB->SetRun(irun);
155   fDButil->UpdateFromCalibDB();
156   fDButil->SetReferenceRun(irun);
157   fDButil->UpdateRefDataFromOCDB();
158   fCalibDB->CreateGUITree("calPads.root");
159   fDButil->CreateGUIRefTree("calPadsRef.root");
160
161   //
162   AliSplineFit *fitVdrift=0x0;
163   Int_t startTimeGRP=0, stopTimeGRP=0;
164   if (fCalibDB->GetGRP(irun)){
165     startTimeGRP = AliTPCcalibDB::GetGRP(irun)->GetTimeStart();
166     stopTimeGRP  = AliTPCcalibDB::GetGRP(irun)->GetTimeEnd();
167   }
168   if (startTime==0){
169     startTime=startTimeGRP;
170     endTime=stopTimeGRP;
171   }
172   AliTPCSensorTempArray * tempArray = fCalibDB->GetTemperatureSensor(irun);
173   AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
174   AliDCSSensorArray* goofieArray = fCalibDB->GetGoofieSensors(irun);
175   //
176   Int_t dtime = TMath::Max((endTime-startTime)/20,10*60);
177   //
178   //Goofie statistical data
179   //
180   TVectorD vecEntries, vecMean, vecMedian,vecRMS;
181   fDButil->ProcessGoofie(vecEntries ,vecMedian, vecMean, vecRMS);
182   //
183   //CE data processing - see ProcessCEdata function for description of the results
184   //
185   TVectorD fitResultsA, fitResultsC;
186   Int_t nmaskedCE;
187   Double_t chi2ACE=0,chi2CCE=0;
188   //     fDButil->ProcessCEdata("(sector<36)++gx++gy++lx++lx**2",fitResultsA,fitResultsC,nmaskedCE);
189   fDButil->ProcessCEdata("(sector<36)++gx++gy++(lx-134)++(sector<36)*(lx-134)++(ly/lx)^2",fitResultsA,fitResultsC,nmaskedCE,chi2ACE,chi2CCE);
190
191   TVectorD fitCEResultsA(7), fitCEResultsC(7);
192   Int_t    noutCE;
193   Double_t chi2CEA=0,chi2CEC=0;
194   AliTPCCalPad *time0 = fDButil->CreatePadTime0CE(fitCEResultsA, fitCEResultsC, noutCE, chi2CEA, chi2CEC);
195   delete time0;
196   //  
197   //
198   TVectorD vecTEntries, vecTMean, vecTRMS, vecTMedian, vecQEntries, vecQMean, vecQRMS, vecQMedian;
199   Float_t driftTimeA, driftTimeC;
200   fDButil->ProcessCEgraphs(vecTEntries, vecTMean, vecTRMS, vecTMedian,
201                           vecQEntries, vecQMean, vecQRMS, vecQMedian,
202                           driftTimeA, driftTimeC );
203   //
204   //
205   //
206   //drift velocity using tracks
207   //
208   //     fitVdrift=fCalibDB->GetVdriftSplineFit("ALISPLINEFIT_MEAN_VDRIFT_COSMICS_ALL",irun);
209   fitVdrift=fCalibDB->CreateVdriftSplineFit("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL",irun);
210   //noise data Processing - see ProcessNoiseData function for description of the results
211   TVectorD vNoiseMean, vNoiseMeanSenRegions, vNoiseRMS, vNoiseRMSSenRegions;
212   Int_t nonMaskedZero=0, nNaN=0;
213   fDButil->ProcessNoiseData(vNoiseMean, vNoiseMeanSenRegions, vNoiseRMS, vNoiseRMSSenRegions, nonMaskedZero, nNaN);
214   //
215   // comparisons
216   //
217   TVectorF pedestalDeviations;
218   TVectorF noiseDeviations;
219   TVectorF pulserQdeviations;
220   Float_t varQMean;
221   Int_t npadsOutOneTB;
222   Int_t npadsOffAdd;
223   fDButil->ProcessPedestalVariations(pedestalDeviations);
224   fDButil->ProcessNoiseVariations(noiseDeviations);
225   fDButil->ProcessPulserVariations(pulserQdeviations,varQMean,npadsOutOneTB,npadsOffAdd);
226   //
227   //L3 data 
228   //
229   Float_t bz=AliTPCcalibDB::GetBz(irun);
230   Char_t  l3pol=AliTPCcalibDB::GetL3Polarity(irun);
231   //
232   //QA data processing
233   //
234   TVectorD vQaOcc;
235   TVectorD vQaQtot;
236   TVectorD vQaQmax;
237   fDButil->ProcessQAData(vQaOcc, vQaQtot, vQaQmax);
238   //
239   //calibration Pulser data processing
240   //
241   Int_t nOffChannels=0;
242   TVectorD vTimePulser;
243   nOffChannels=fDButil->GetNPulserOutliers();
244   fDButil->ProcessPulser(vTimePulser);
245   //
246   //ALTRO data
247   //
248   Int_t nMasked=0;
249   fDButil->ProcessALTROConfig(nMasked);
250   //
251   //Calib RAW data
252   //
253   Int_t nFailL1=-1;
254   if (fCalibDB->GetCalibRaw()) nFailL1=fCalibDB->GetCalibRaw()->GetNFailL1Phase();
255   //
256   //production information
257   //
258   Int_t nalien=0,nRawAlien=0,nlocal=0,nRawLocal=0;
259   //run type
260   TObjString runType(AliTPCcalibDB::GetRunType(irun).Data());
261   //
262   //
263   //
264   
265   for (Int_t itime=startTime; itime<endTime; itime+=dtime){
266     //
267     TTimeStamp tstamp(itime);
268     Float_t valuePressure  = fCalibDB->GetPressure(tstamp,irun,0);
269     Float_t valuePressure2 = fCalibDB->GetPressure(tstamp,irun,1);
270     Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative((UInt_t)itime,irun,0);
271     Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative((UInt_t)itime,irun,1);
272     //temperature fits
273     TLinearFitter * fitter = 0;
274     TVectorD vecTemp[10];
275     for (Int_t itype=0; itype<5; itype++)
276       for (Int_t iside=0; iside<2; iside++){
277         fitter= tempMap->GetLinearFitter(itype,iside,tstamp);     
278         if (!fitter) continue;
279         fitter->Eval();
280         fitter->GetParameters(vecTemp[itype+iside*5]);
281         delete fitter;
282       }
283     //
284     //measured skirt temperatures
285     //
286     TVectorD vecSkirtTempA(18);
287     TVectorD vecSkirtTempC(18);
288     Int_t nsenTemp=tempArray->NumSensors();
289     for (Int_t isenTemp=0;isenTemp<nsenTemp;++isenTemp){
290       AliTPCSensorTemp *senTemp=(AliTPCSensorTemp*)tempArray->GetSensorNum(isenTemp);
291       if (senTemp->GetType()!=3) continue;
292       if (TMath::Sqrt(senTemp->GetX()*senTemp->GetX()+senTemp->GetY()*senTemp->GetY())<100) continue; //only skirt, outer FC vessel
293       Double_t val=senTemp->GetValue(tstamp);
294       if (senTemp->GetSide()==0)
295         vecSkirtTempA[senTemp->GetSector()]=val;
296       else
297         vecSkirtTempC[senTemp->GetSector()]=val;
298     }
299     //
300     //goofie data
301     //
302     TVectorD vecGoofie; 
303     if (goofieArray){
304       vecGoofie.ResizeTo(goofieArray->NumSensors());
305       for (Int_t isensor=0; isensor<goofieArray->NumSensors();isensor++){
306         AliDCSSensor *gsensor = goofieArray->GetSensor(isensor);
307         if (gsensor){
308           vecGoofie[isensor] = gsensor->GetValue(tstamp);
309         }
310       }
311     } else {
312       vecGoofie.ResizeTo(19);
313     }
314     //
315     TVectorD voltagesIROC(36);
316     TVectorD voltagesOROC(36);
317     for(Int_t j=1; j<36; j++) voltagesIROC[j-1] = fCalibDB->GetChamberHighVoltage(irun, j,itime);
318     for(Int_t j=36; j<72; j++) voltagesOROC[j-36] = fCalibDB->GetChamberHighVoltage(irun, j,itime);
319     Double_t voltIROC = TMath::Median(36, voltagesIROC.GetMatrixArray());
320     Double_t voltOROC = TMath::Median(36, voltagesOROC.GetMatrixArray());
321     //
322     Float_t  coverIA=AliTPCcalibDB::GetCoverVoltage(irun,0,itime);
323     Float_t  coverIC=AliTPCcalibDB::GetCoverVoltage(irun,18,itime);
324     Float_t  coverOA=AliTPCcalibDB::GetCoverVoltage(irun,36,itime);
325     Float_t  coverOC=AliTPCcalibDB::GetCoverVoltage(irun,54,itime);
326     Float_t  skirtA=AliTPCcalibDB::GetSkirtVoltage(irun,0,itime);
327     Float_t  skirtC=AliTPCcalibDB::GetSkirtVoltage(irun,18,itime);
328     Float_t  ggOffA=AliTPCcalibDB::GetGGoffsetVoltage(irun,0,itime);
329     Float_t  ggOffC=AliTPCcalibDB::GetGGoffsetVoltage(irun,18,itime);
330     //drift velocity
331     Float_t dvCorr=-5;
332     if (fitVdrift) dvCorr=fitVdrift->Eval(itime);
333     
334     //tempMap->GetLinearFitter(0,0,itime);
335     (*fPcstream)<<"dcs"<<
336       "run="<<irun<<
337       "time="<<itime<<
338       "startTimeGRP="<<startTimeGRP<<
339       "stopTimeGRP="<<stopTimeGRP<<
340       //run type
341       "runType.="<<&runType<<
342       // voltage setting
343       "VIROC.="<<&voltagesIROC<<
344       "VOROC.="<<&voltagesOROC<<
345       "medianVIROC="<<voltIROC<<
346       "medianVOROC="<<voltOROC<<
347       "coverIA=" << coverIA <<
348       "coverIC=" << coverIC <<
349       "coverOA=" << coverOA <<
350       "coverOC=" << coverOC <<
351       "skirtA=" << skirtA <<
352       "skirtC=" << skirtC <<
353       "ggOffA=" << ggOffA <<
354       "ggOffC=" << ggOffC <<
355       //
356       "ptrel0="<<ptrelative0<<  // deltaTP/TP  - A side
357       "ptrel1="<<ptrelative1<<  // deltaTP/TPC - C side
358       "goofie.="<<&vecGoofie<<
359       "goofieE.="<<&vecEntries<<
360       "goofieMean.="<<&vecMean<<
361       "goofieMedian.="<<&vecMedian<<
362       "goofieRMS.="<<&vecRMS<<
363       //
364       "press="<<valuePressure<<
365       "press2="<<valuePressure2<<
366       "temp00.="<<&vecTemp[0]<<
367       "temp10.="<<&vecTemp[1]<<
368       "temp20.="<<&vecTemp[2]<<
369       "temp30.="<<&vecTemp[3]<<
370       "temp40.="<<&vecTemp[4]<<
371       "temp01.="<<&vecTemp[5]<<
372       "temp11.="<<&vecTemp[6]<<
373       "temp21.="<<&vecTemp[7]<<
374       "temp31.="<<&vecTemp[8]<<
375       "temp41.="<<&vecTemp[9]<<
376       "tempSkirtA.="<<&vecSkirtTempA<<
377       "tempSkirtC.="<<&vecSkirtTempC;
378
379     ProcessDrift(irun, itime);
380     ProcessDriftCE(irun,itime);
381     ProcessDriftAll(irun,itime);
382     //    ProcessKryptonTime(irun,itime);
383     ProcessCTP(irun,itime);
384     ProcessAlign(irun,itime);
385     ProcessGain(irun,itime);
386     ProcessDriftCERef();
387     ProcessPulserRef();
388
389
390     (*fPcstream)<<"dcs"<<       
391       //noise data
392       "meanNoise.="<<&vNoiseMean<<
393       "meanNoiseSen.="<<&vNoiseMeanSenRegions<<
394       "rmsNoise.="<<&vNoiseRMS<<
395       "rmsNoiseSen.="<<&vNoiseRMSSenRegions<<
396       "zeroNoise="<<nonMaskedZero<<
397       "nNaN="<<nNaN<<  
398       //QA data
399       "occQA.="  << &vQaOcc  <<
400       "qQA.="    << &vQaQtot <<
401       "qmaxQA.=" << &vQaQmax <<
402       //pulser data
403       "timePulser.=" << &vTimePulser <<
404       "nOffPulser="<<nOffChannels<<
405       //altro data
406       "nMasked="<<nMasked<<
407       //ce data -Jens version
408       "CEfitA.="<<&fitResultsA<<
409       "CEfitC.="<<&fitResultsC<<
410       "nmaskedCE="<<nmaskedCE<<
411       "chi2ACE="<<chi2ACE<<
412       "chi2CCE="<<chi2CCE<<
413       //ce data new - MI version
414       "CEfitAMI.="<<&fitCEResultsA<<
415       "CEfitCMI.="<<&fitCEResultsC<<
416       "chi2CEA="<<chi2CEA<<
417       "chi2CEC="<<chi2CEC<<
418       //
419       //ce graph data
420       "CEgrTEntries.="<<&vecTEntries<<
421       "CEgrTMean.="<<&vecTMean<<
422       "CEgrTRMS.="<<&vecTRMS<<
423       "CEgrTMedian.="<<&vecTMedian<<
424       "CEgrQEntries.="<<&vecQEntries<<
425       "CEgrQMean.="<<&vecQMean<<
426       "CEgrQRMS.="<<&vecQRMS<<
427       "CEgrQMedian.="<<&vecQMedian<<
428       "CEgrDriftA="<<driftTimeA<<
429       "CEgrDriftC="<<driftTimeC<<
430       //calib raw data
431       "nFailL1="<<nFailL1<<
432       // b field
433       "Bz="<< bz <<
434       "L3polarity="<<l3pol<<
435       // production information
436       "nalien="<<nalien<<
437       "nRawAlien="<<nRawAlien<<
438       "nlocal="<<nlocal<<
439       "nRawLocal="<<nRawLocal<<
440       //comparisons with ref data
441       "pedestalDeviations.="<<&pedestalDeviations<<
442       "noiseDeviations.="<<&noiseDeviations<<
443       "pulserQdeviations.="<<&pulserQdeviations<<
444       //         "pulserVarQMean="<<varQMean<<
445       "pulserNpadsOutOneTB="<<npadsOutOneTB<<
446       "pulserNpadsOffAdd="<<npadsOffAdd<<
447       "driftCorrCosmAll="<<dvCorr<<
448       "\n";
449   }//end run loop
450 }
451
452
453
454
455
456
457 void AliTPCcalibSummary::ProcessDrift(Int_t run, Int_t timeStamp){
458   //
459   // dump drift calibration data to the tree
460   //
461   TObjArray *array =fCalibDB->GetTimeVdriftSplineRun(run);
462   TGraphErrors *laserA[3]={0,0,0};
463   TGraphErrors *laserC[3]={0,0,0};
464   TGraphErrors *cosmicAll=0;
465   static Double_t     vlaserA[3]={0,0,0};
466   static Double_t     vlaserC[3]={0,0,0};
467   static Double_t     vcosmicAll=0;
468   static Double_t     vdrift1=0;
469   vdrift1=fCalibDB->GetVDriftCorrectionTime(timeStamp,run,0,1);
470
471   if (array){
472     laserA[0]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DELAY_LASER_ALL_A");
473     laserA[1]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_A");
474     laserA[2]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
475     laserC[0]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DELAY_LASER_ALL_C");
476     laserC[1]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_C");
477     laserC[2]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
478     cosmicAll =(TGraphErrors*)array->FindObject("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL");
479   }
480   if (laserA[0]) vlaserA[0]= AliTPCcalibDButil::EvalGraphConst(laserA[0],timeStamp);
481   if (laserA[1]) vlaserA[1]= AliTPCcalibDButil::EvalGraphConst(laserA[1],timeStamp);
482   if (laserA[2]) vlaserA[2]= AliTPCcalibDButil::EvalGraphConst(laserA[2],timeStamp);
483   if (laserC[0]) vlaserC[0]= AliTPCcalibDButil::EvalGraphConst(laserC[0],timeStamp);
484   if (laserC[1]) vlaserC[1]= AliTPCcalibDButil::EvalGraphConst(laserC[1],timeStamp);
485   if (laserC[2]) vlaserC[2]= AliTPCcalibDButil::EvalGraphConst(laserC[2],timeStamp);
486   if (cosmicAll) vcosmicAll= AliTPCcalibDButil::EvalGraphConst(cosmicAll,timeStamp); 
487   (*fPcstream)<<"dcs"<<
488     "vlaserA0="<<vlaserA[0]<<
489     "vlaserA1="<<vlaserA[1]<<
490     "vlaserA2="<<vlaserA[2]<<
491     "vlaserC0="<<vlaserC[0]<<
492     "vlaserC1="<<vlaserC[1]<<
493     "vlaserC2="<<vlaserC[2]<<
494     "vcosmicAll="<<vcosmicAll<<
495     "vdrift1="<<vdrift1;
496
497   //
498   // define distance to measurement
499   //
500   static Double_t dlaserA=0; 
501   static Double_t dlaserC=0; 
502   static Double_t dcosmic=0; 
503   static Double_t slaserA=0; 
504   static Double_t slaserC=0; 
505   static Double_t scosmic=0; 
506   static Double_t  vclaserA[3]={0,0,0};
507   static Double_t  vclaserC[3]={0,0,0};
508   static Double_t  vccosmicAll=0;
509   for (Int_t i=0;i<3;i++){
510     if (laserA[i]) AliTPCcalibDButil::GetNearest(laserA[i],timeStamp,dlaserA,vclaserA[i]);
511     if (laserC[i]) AliTPCcalibDButil::GetNearest(laserC[i],timeStamp,dlaserC,vclaserC[i]);
512   }  
513   if (cosmicAll) AliTPCcalibDButil::GetNearest(cosmicAll,timeStamp,dcosmic,vccosmicAll);
514   (*fPcstream)<<"dcs"<<
515     "vclaserA0="<<vclaserA[0]<<
516     "vclaserA1="<<vclaserA[1]<<
517     "vclaserA2="<<vclaserA[2]<<
518     "vclaserC0="<<vclaserC[0]<<
519     "vclaserC1="<<vclaserC[1]<<
520     "vclaserC2="<<vclaserC[2]<<
521     "vccosmicAll="<<vccosmicAll<<
522     "dlaserA="<<dlaserA<<
523     "dlaserC="<<dlaserC<<
524     "dcosmic="<<dcosmic<<
525     "slaserA="<<slaserA<<
526     "slaserC="<<slaserC<<
527     "scosmic="<<scosmic;
528 }
529
530 void AliTPCcalibSummary::ProcessDriftCE(Int_t run,Int_t timeStamp){
531   //
532   // dump drift calibration data CE
533   //
534   TObjArray *arrT=fCalibDB->GetCErocTtime();
535   AliTPCParam *param=fCalibDB->GetParameters();
536   static TVectorD tdriftCE(74);
537   static TVectorD tndriftCE(74);
538   static TVectorD vdriftCE(74);
539   static TVectorD tcdriftCE(74);
540   static TVectorD tddriftCE(74);
541   static Double_t ltime0A;
542   static Double_t ltime0C;
543   //
544   //
545   //
546   ltime0A  = fDButil->GetLaserTime0(run,timeStamp,36000,0);
547   ltime0C  = fDButil->GetLaserTime0(run,timeStamp,36000,1);
548   //
549   for (Int_t i=0; i<arrT->GetEntries();i++){
550     tdriftCE[i]=0;
551     vdriftCE[i]=0;
552     TGraph *graph = (TGraph*)arrT->At(i);
553     if (!graph) continue;
554     tdriftCE[i]=AliTPCcalibDButil::EvalGraphConst(graph,timeStamp);
555     Double_t deltaT,gry;
556     AliTPCcalibDButil::GetNearest(graph,timeStamp,deltaT,gry);
557     tndriftCE[i]=graph->GetN();
558     tcdriftCE[i]=gry;          
559     tddriftCE[i]=deltaT;               
560     if (i%36<18){
561       vdriftCE[i] =(param->GetZLength(i)/(tdriftCE[i]*param->GetTSample()*(1.-ltime0A)-param->GetL1Delay()))/param->GetDriftV();
562     }else{
563       vdriftCE[i] =(param->GetZLength(i)/(tdriftCE[i]*param->GetTSample()*(1.-ltime0A)-param->GetL1Delay()))/param->GetDriftV();
564     }
565   }
566   // export values
567   (*fPcstream)<<"dcs"<<  
568     "tdriftCE.="<<&tdriftCE<<    // arrival time 
569     "vdriftCE.="<<&vdriftCE<<    // derived drift velocity per chamber 
570     "tndriftCE.="<<&tndriftCE<<  // number of points (chambers)
571     "tcdriftCE.="<<&tcdriftCE<<  // constant evaluation - nearest point used
572     "tddriftCE.="<<&tddriftCE<<  // distance to closest measuement
573     "ltime0A="<<ltime0A<<        // laser offset expected in reconstruction
574     "ltime0C="<<ltime0C;         // laser offset expected in reconstruction 
575    }
576
577
578 void AliTPCcalibSummary::ProcessDriftAll(Int_t run,Int_t timeStamp){
579   //
580   // dump drift calibration data  all calibrations form DB util
581   // test of utils
582   static Double_t vdriftCEA=0, vdriftCEC=0, vdriftCEM=0;
583   static Double_t vdriftLTA=0, vdriftLTC=0, vdriftLTM=0;
584   static Double_t vdriftITS=0;
585   static Double_t vdriftP=0;
586   static Double_t dcea=0, dcec=0, dcem=0,  dla=0,dlc=0,dlm=0,dp=0;
587   static Double_t dits=0;
588   static Double_t ltime0A;
589   static Double_t ltime0C;
590   static Double_t ctime0;
591   static Double_t vdrift1=0;
592   vdrift1= fCalibDB->GetVDriftCorrectionTime(timeStamp,run,0,1);
593   vdriftP = fDButil->GetVDriftTPC(dp, run, timeStamp, 86400, 3600,0);
594   ctime0 = AliTPCcalibDButil::GetTriggerOffsetTPC(run,timeStamp, 36000, 3600,0);
595   //
596   vdriftCEA= fDButil->GetVDriftTPCCE(dcea,run,timeStamp,36000,0);
597   vdriftCEC= fDButil->GetVDriftTPCCE(dcec,run,timeStamp,36000,1);
598   vdriftCEM= fDButil->GetVDriftTPCCE(dcem,run,timeStamp,36000,2);
599   //
600   vdriftLTA= fDButil->GetVDriftTPCLaserTracks(dla,run,timeStamp,36000,0);
601   vdriftLTC= fDButil->GetVDriftTPCLaserTracks(dlc,run,timeStamp,36000,1);
602   vdriftLTM= fDButil->GetVDriftTPCLaserTracks(dlm,run,timeStamp,36000,2);
603   //
604   vdriftITS= fDButil->GetVDriftTPCITS(dits, run,timeStamp);
605   //
606   ltime0A  = fDButil->GetLaserTime0(run,timeStamp,36000,0);
607   ltime0C  = fDButil->GetLaserTime0(run,timeStamp,36000,1);
608
609   (*fPcstream)<<"dcs"<<  
610     //
611     "vdriftCEA="<<vdriftCEA<<   // drift velocity CE
612     "vdriftCEC="<<vdriftCEC<<
613     "vdriftCEM="<<vdriftCEM<<
614     "dcea="<<dcea<<
615     "dcec="<<dcec<<
616     "dcem="<<dcem<<
617     "vdriftLTA="<<vdriftLTA<<   // drift velocity laser tracks
618     "vdriftLTC="<<vdriftLTC<<
619     "vdriftLTM="<<vdriftLTM<<
620     "dla="<<dla<<
621     "dlc="<<dlc<<
622     "dlm="<<dlm<<
623     //
624     //
625     "vdriftITS="<<vdriftITS<<
626     "dits="<<dits<<
627     "ctime0="<<ctime0<<
628     "vdriftP="<<vdriftP<<       // drift velocity comsic 
629     "dp="<<dp<<
630     "vdrift1="<<vdrift1;        // combined drift velocity
631
632 }
633
634
635
636 void AliTPCcalibSummary::ProcessKryptonTime(Int_t run, Int_t timeStamp){
637   //
638   // Dumping  krypton calibration results
639   //
640   static TObjArray * krArray=0;
641   if (!krArray) {
642     AliCDBEntry* entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeGainKrypton", run);
643     if (entry){
644       krArray = (TObjArray*)entry->GetObject();
645     }
646   }
647   static TVectorD krMean(74);
648   static TVectorD krErr(74);
649   static TVectorD krDist(74);
650   TGraphErrors *gr=0;
651   Double_t deltaT=0,gry=0;
652   if (krArray){
653     for (Int_t isec=0; isec<72; isec++){
654       krMean[isec]=0;
655       krDist[isec]=0;
656       krErr[isec]=0;
657       gr=(TGraphErrors*)krArray->At(isec);
658       if (gr) {
659         krMean[isec]=AliTPCcalibDButil::EvalGraphConst(gr,timeStamp);
660         AliTPCcalibDButil::GetNearest(gr, timeStamp,deltaT,gry);
661         krDist[isec]=deltaT;
662       }     
663       if (72+isec<krArray->GetEntries()) {
664         gr=(TGraphErrors*)krArray->At(72+isec);
665         if (gr) krErr[isec]=AliTPCcalibDButil::EvalGraphConst(gr,timeStamp);
666       }
667     }
668     krMean[72]= TMath::Median(36,krMean.GetMatrixArray());
669     krMean[73]= TMath::Median(36,&(krMean.GetMatrixArray()[36]));
670     krErr[72]= TMath::Median(36,krErr.GetMatrixArray());
671     krErr[73]= TMath::Median(36,&(krErr.GetMatrixArray()[36]));
672   }
673   (*fPcstream)<<"dcs"<<
674     "krMean.="<<&krMean<<
675     "krErr.="<<&krErr<<
676     "krDist.="<<&krDist;
677 }
678
679 void AliTPCcalibSummary::ProcessCTP(Int_t irun, Int_t timeStamp){
680   //
681   // 
682   //
683   static TClonesArray *pcarray = new TClonesArray("AliCTPInputTimeParams",1);
684   static AliCTPTimeParams *ctpParams =0;
685   ctpParams = fCalibDB->GetCTPTimeParams(); //
686   const TObjArray        *arr = ctpParams->GetInputTimeParams();
687   pcarray->ExpandCreateFast(TMath::Max(arr->GetEntries(),1));
688   for (Int_t i=0; i<arr->GetEntries(); i++){
689     new ((*pcarray)[i]) AliCTPInputTimeParams(*((AliCTPInputTimeParams*)arr->At(i)));
690   }
691   (*fPcstream)<<"ctp"<<
692     "run="<<irun<<
693     "time="<<timeStamp<<
694     "ctpP.="<<ctpParams<<
695     "ctpArr="<<pcarray<<
696     "\n";
697   (*fPcstream)<<"dcs"<<
698     "ctpP.="<<ctpParams<<
699     "ctpArr="<<pcarray;
700 }
701
702 void AliTPCcalibSummary::ProcessGain(Int_t irun, Int_t timeStamp){
703   //
704   // Dump gain related information to the tree
705   //
706   static Float_t  gainCosmic = 0;
707   static Float_t  gainMIP = 0;
708   static Float_t  attachMIP = 0;
709   static Double_t dMIP=0; 
710   Double_t dummy=0;
711   TObjArray * gainSplines = fCalibDB->GetTimeGainSplinesRun(irun);
712   if (gainSplines) {
713     TGraphErrors * graphMIP = (TGraphErrors *) gainSplines->FindObject("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");
714     TGraphErrors * graphCosmic = (TGraphErrors *) gainSplines->FindObject("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL");
715     TGraphErrors * graphAttach = (TGraphErrors *) gainSplines->FindObject("TGRAPHERRORS_MEAN_ATTACHMENT_BEAM_ALL");
716     if (graphMIP) gainMIP = AliTPCcalibDButil::EvalGraphConst(graphMIP,timeStamp);
717     if (graphCosmic) gainCosmic = AliTPCcalibDButil::EvalGraphConst(graphCosmic,timeStamp);
718     if (graphAttach) attachMIP = AliTPCcalibDButil::EvalGraphConst(graphAttach,timeStamp);
719     if (graphMIP)  AliTPCcalibDButil::GetNearest(graphMIP, timeStamp, dMIP,dummy);
720   }
721   // time dependence of gain
722   (*fPcstream)<<"dcs"<<
723     "gainMIP="<<gainMIP<<
724     "attachMIP="<<attachMIP<<
725     "dMIP="<<dMIP<<
726     "gainCosmic="<<gainCosmic;     
727 }
728
729
730 void AliTPCcalibSummary::ProcessAlign(Int_t run, Int_t timeStamp){
731   //
732   // Proccess alignment 
733   //
734   TString   grnames[12]={"ALIGN_ITS", "ALIGN_ITSP", "ALIGN_ITSM", "ALIGN_ITSB",
735                          "ALIGN_TRD", "ALIGN_TRDP", "ALIGN_TRDM","ALIGN_TRDB",
736                          "ALIGN_TOF", "ALIGN_TOFP", "ALIGN_TOFM","ALIGN_TOFB"};
737   TString   grpar[9]={"DELTAPSI", "DELTATHETA", "DELTAPHI",
738                       "DELTAX", "DELTAY", "DELTAZ",
739                       "DRIFTVD", "T0", "VDGY"};
740   static Double_t values[12*9];
741   static Double_t errs[12*9];
742
743   TObjArray *array =fCalibDB->GetTimeVdriftSplineRun(run);
744   TGraphErrors *gr=0;
745   for (Int_t idet=0; idet<12; idet++){
746     for (Int_t ipar=0; ipar<9; ipar++){
747       TString grName=grnames[idet];
748       grName+="_TPC_";
749       grName+=grpar[ipar];
750       if (array){
751         gr = (TGraphErrors*)array->FindObject(grName.Data());
752       }
753       values[9*idet+ipar]=0;
754       errs[9*idet+ipar]=0;
755       if (gr) values[9*idet+ipar]=AliTPCcalibDButil::EvalGraphConst(gr,timeStamp);
756       (*fPcstream)<<"dcs"<<
757         Form("%s=",grName.Data())<<values[9*idet+ipar];
758       (*fPcstream)<<"align"<<
759         Form("%s=",grName.Data())<<values[9*idet+ipar];
760     }
761   }
762   (*fPcstream)<<"align"<<
763     "time="<<timeStamp<<
764     "run="<<run<<
765     "\n";
766 }
767
768
769 void AliTPCcalibSummary::ProcessDriftCERef(){
770   //
771   // Get fit of residuals if CE in respect with reference
772   // data
773   //
774   static TVectorD  sec(72);
775   static TVectorD vec0(72);
776   static TVectorD vecLy(72);
777   static TVectorD vecLx(72);
778   static TVectorD vecChi2(72);
779   static TVectorD vecN(72);
780   //
781   static TVectorD vecA0(72);
782   static TVectorD vecALy(72);
783   static TVectorD vecALx(72);
784   static TVectorD vecAChi2(72);
785   //
786   static TVectorD vecASide(4);
787   static TVectorD vecCSide(4);
788   static Bool_t isCalc=kFALSE;
789   
790   TFile f("calPads.root");
791   TFile fref("calPadsRef.root");
792   TTree * tree = (TTree*)f.Get("calPads");
793   TTree * treeRef = (TTree*)fref.Get("calPads");
794   tree->AddFriend(treeRef,"R");
795   tree->SetAlias("inCE","((CEQmean.fElements>35)&&abs(CETmean.fElements)<1.5&&abs(CETrms.fElements/1.2-1)<0.2)");  // outlyerTrms
796   tree->SetAlias("inCER","((R.CEQmean.fElements>35)&&abs(R.CETmean.fElements)<1.5&&abs(R.CETrms.fElements/1.2-1)<0.2)");  // outlyerTrms
797   //
798   if (!isCalc){
799     // make fits only once
800     TStatToolkit toolkit;
801     Double_t chi2=0;
802     Int_t    npoints=0;
803     TVectorD param;
804     TMatrixD covar;
805     tree->SetAlias("dt","CETmean.fElements-R.CETmean.fElements");
806     TCut cutAll ="inCE&&inCER&&abs(CETmean.fElements-R.CETmean.fElements)<0.5"; 
807     TString  fstringG="";              // global part
808     fstringG+="ly.fElements++";
809     fstringG+="(lx.fElements-134.)++";  
810     for (Int_t isec=0; isec<72; isec++){
811       TStatToolkit::FitPlane(tree,"0.264*dt", fstringG.Data(),Form("sector==%d",isec)+cutAll, chi2,npoints,param,covar,-1,0, 10000000, kFALSE);
812       if (npoints<3) continue;
813       printf("Sector=%d\n",isec);
814       vec0[isec]=param[0];
815       vecLy[isec]=param[1];
816       vecLx[isec]=param[2];
817       sec[isec]=isec;
818       vecN[isec]=npoints;
819       vecChi2[isec]=TMath::Sqrt(chi2/npoints);
820
821       TStatToolkit::FitPlane(tree,"0.264*CETmean.fElements", fstringG.Data(),Form("sector==%d",isec)+cutAll, chi2,npoints,param,covar,-1,0, 10000000, kFALSE);
822       if (npoints<3) continue;
823       printf("Sector=%d\n",isec);
824       vecA0[isec]=param[0];
825       vecALy[isec]=param[1];
826       vecALx[isec]=param[2];
827       vecAChi2[isec]=TMath::Sqrt(chi2/npoints);
828     }
829     isCalc=kTRUE;
830     //
831     TString  fstringRef="";              // global fit
832     fstringRef+="gx.fElements++";
833     fstringRef+="gy.fElements++";  
834     fstringRef+="lx.fElements++";
835     TStatToolkit::FitPlane(tree,"0.264*dt", fstringG.Data(),"(sector%36)<18"+cutAll, chi2,npoints,vecASide,covar,-1,0, 10000000, kFALSE);
836     TStatToolkit::FitPlane(tree,"0.264*dt", fstringG.Data(),"(sector%36)>=18"+cutAll, chi2,npoints,vecCSide,covar,-1,0, 10000000, kFALSE);
837
838   }
839   (*fPcstream)<<"dcs"<<     // CE information
840     "CETSector.="<<&sec<<    // sector numbers
841     "CETRefA.="<<&vecASide<<   // diff to reference A side
842     "CETRefC.="<<&vecCSide<<   // diff to reference C side    
843     //                      // fit in respect to reference data
844     "CETRef0.="<<&vec0<<    // offset change
845     "CETRefY.="<<&vecLy<<   // slope y change - rad
846     "CETRefX.="<<&vecLx<<   // slope x change - rad
847     "CETRefChi2.="<<&vecChi2<<    // chi2 (rms in cm)
848     "CETRefN.="<<&vecN<<     //number of accepted points
849     //                       // fit in respect per mean per side
850     "CET0.="<<&vecA0<<       // offset change
851     "CETY.="<<&vecALy<<      // slope y change - rad
852     "CETX.="<<&vecALx<<      // slope x change - rad
853     "CETChi2.="<<&vecAChi2;  // chi2 (rms in cm)
854 }
855
856 void AliTPCcalibSummary::ProcessPulserRef(){
857   //
858   // Get fit of residuals if Pulser in respect with reference
859   // data
860   //
861   static TVectorD  sec(72);
862   static TVectorD vec0(72);
863   static TVectorD vecLy(72);
864   static TVectorD vecLx(72);
865   static TVectorD vecChi2(72);
866   static TVectorD vecN(72);
867   //
868   static TVectorD vecA0(72);
869   static TVectorD vecALy(72);
870   static TVectorD vecALx(72);
871   static TVectorD vecAChi2(72);
872   static Bool_t isCalc=kFALSE;
873   
874   TFile f("calPads.root");
875   TFile fref("calPadsRef.root");
876   TTree * tree = (TTree*)f.Get("calPads");
877   TTree * treeRef = (TTree*)fref.Get("calPads");
878   tree->AddFriend(treeRef,"R");
879   
880   tree->SetAlias("inPulser","(abs(PulserTmean.fElements-PulserTmean_Median)<1.5&&abs(PulserTrms.fElements-PulserTrms_Median)<0.2)");  // outlyerTrms
881   tree->SetAlias("inPulserR","(abs(R.PulserTmean.fElements-R.PulserTmean_Median)<1.5&&abs(R.PulserTrms.fElements-R.PulserTrms_Median)<0.2)");  // outlyerTrms
882   //
883   if (!isCalc){
884     // make fits only once
885     TStatToolkit toolkit;
886     Double_t chi2=0;
887     Int_t    npoints=0;
888     TVectorD param;
889     TMatrixD covar;
890     tree->SetAlias("dt","PulserTmean.fElements-R.PulserTmean.fElements");
891     TCut cutAll ="inPulser&&inPulserR"; 
892     TString  fstringG="";              // global part
893     fstringG+="ly.fElements++";
894     fstringG+="(lx.fElements-134.)++";  
895     for (Int_t isec=0; isec<72; isec++){
896       TStatToolkit::FitPlane(tree,"dt", fstringG.Data(),Form("sector==%d",isec)+cutAll, chi2,npoints,param,covar,-1,0, 10000000, kFALSE);
897       if (npoints<3) continue;
898       printf("Setor=%d\n",isec);
899       vec0[isec]=param[0];
900       vecLy[isec]=param[1];
901       vecLx[isec]=param[2];
902       sec[isec]=isec;
903       vecN[isec]=npoints;
904
905       TStatToolkit::FitPlane(tree,"PulserTmean.fElements", fstringG.Data(),Form("sector==%d",isec)+cutAll, chi2,npoints,param,covar,-1,0, 10000000, kFALSE);
906       if (npoints<3) continue;
907       printf("Setor=%d\n",isec);
908       vecA0[isec]=param[0];
909       vecALy[isec]=param[1];
910       vecALx[isec]=param[2];
911       vecAChi2[isec]=TMath::Sqrt(chi2/npoints);
912     }
913
914     isCalc=kTRUE;
915   }
916   (*fPcstream)<<"dcs"<<     // Pulser information
917     "PulserTSector.="<<&sec<<    // sector numbers
918     //                      // fit in respect to reference
919     "PulserTRef0.="<<&vec0<<    // offset change
920     "PulserTRefY.="<<&vecLy<<   // slope y change - rad
921     "PulserTRefX.="<<&vecLx<<   // slope x change - rad
922     "PulserTRefChi2.="<<&vecChi2<<    // chi2 (rms in cm)
923     "PulserTRefN.="<<&vecN<<     //number of accepted points
924     //                       // fit in respect per mean per side
925     "PulserT0.="<<&vecA0<<       // offset change
926     "PulserTY.="<<&vecALy<<      // slope y change - rad
927     "PulserTX.="<<&vecALx<<      // slope x change - rad
928     "PulserTChi2.="<<&vecAChi2;  // chi2 (rms in cm)
929 }
930
931
932
933
934
935 // TCanvas * DrawCEDiff(TTree * tree){
936   
937 //   TCanvas *canvasIO = new TCanvas("canvasCEIO","canvasCEIO");
938 //   canvasIO->Divide(6,6);
939 //   for (Int_t isec=0; isec<36; isec++){
940 //     canvasIO->cd(isec+1);
941 //     dcs->Draw(Form("CET0.fElements[%d]-CET0.fElements[%d]",isec+36,isec),Form("abs(CETRef0.fElements[%d])<0.3",isec),"");
942 //     printf("%d\t%f\t%f\n",isec,dcs->GetHistogram()->GetMean(),dcs->GetHistogram()->GetRMS());
943 //   }
944
945 // }