]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/CalibMacros/CalibEnv.C
AliTPCcalibDB.cxx - Calculate the distance to the closest measurement
[u/mrichter/AliRoot.git] / TPC / CalibMacros / CalibEnv.C
1 /*
2 //Make a tree dump of TPC calibration:
3
4 gSystem->AddIncludePath("-I$ALICE_ROOT/TPC");
5 .L $ALICE_ROOT/TPC/CalibMacros/CalibEnv.C+
6
7 CalibEnv("run.list");
8 TFile f("dcsTime.root")
9 */
10 #include "TMVA/TSpline1.h"
11 #include "TMVA/TSpline2.h"
12 #include <iostream>
13 #include <fstream>
14 #include <stdio.h>
15 #include <AliCDBManager.h>
16 #include <AliCDBEntry.h>
17 #include <AliLog.h>
18 #include <AliMagF.h>
19 #include "AliTPCcalibDB.h"
20 #include "AliTPCcalibDButil.h"
21 #include "AliTPCAltroMapping.h"
22 #include "AliTPCExB.h"
23 #include "AliTPCCalROC.h"
24 #include "AliTPCCalPad.h"
25 #include "AliTPCSensorTempArray.h"
26 #include "AliGRPObject.h"
27 #include "AliTPCTransform.h"
28 #include "TFile.h"
29 #include "TKey.h"
30 #include "TObjArray.h"
31 #include "TObjString.h"
32 #include "TString.h"
33 #include "AliTPCCalPad.h"
34 #include "AliTPCROC.h"
35 #include "AliTPCParam.h"
36 #include "AliTPCCalibPulser.h"
37 #include "AliTPCCalibPedestal.h"
38 #include "AliTPCCalibCE.h"
39 #include "AliTPCExBFirst.h"
40 #include "TTreeStream.h"
41 #include "AliTPCTempMap.h"
42 #include "TVectorD.h"
43 #include "TMatrixD.h"
44 #include "AliTPCCalibRaw.h"
45 #include "AliSplineFit.h"
46 #include "TGraphErrors.h"
47
48 //
49 AliTPCcalibDB     *calibDB=0;
50 AliTPCcalibDButil *dbutil =0;
51 TTree * dcsTree=0;
52 TString refFile="dummy.root"; 
53 TTreeSRedirector *pcstream =0;
54 //
55 //
56 void ProcessRun(Int_t irun, Int_t startTime, Int_t endTime);
57 void GetProductionInfo(Int_t run, Int_t &nalien, Int_t &nRawAlien, Int_t &nlocal, Int_t &nRawLocal);
58 void ProcessDrift(Int_t run, Int_t timeStamp);
59 void ProcessDriftCE(Int_t run, Int_t timeStamp);
60 void ProcessDriftAll(Int_t run, Int_t timeStamp);
61 void ProcessKryptonTime(Int_t run, Int_t timeStamp);
62 void CalibEnv(const char * runList, Int_t first=1, Int_t last=-1){
63   //
64   // runList - listOfRuns to process
65   // first   - first run to process
66   // last    - last  to process
67   // 
68   refFile=gSystem->Getenv("REF_DATA_FILE");
69   calibDB = AliTPCcalibDB::Instance();
70   dbutil= new AliTPCcalibDButil;
71   //
72   // make list of runs
73   //
74   ifstream in;
75   in.open(runList);
76   Int_t irun=0;
77   TArrayI runArray(100000);
78   Int_t indexes[100000];
79   Int_t nruns=0;
80   while(in.good()) {
81     in >> irun;
82     if (in.eof()) break;
83     if (irun<first) continue;  // process only subset of list
84     if (last>0 && irun>=last) continue;  // process only subset of list
85     runArray[nruns]=irun;
86     nruns++;
87   }
88   TMath::Sort(nruns, runArray.fArray, indexes,kFALSE);
89   pcstream = new TTreeSRedirector("dcsTime.root");
90   dbutil->SetRefFile(refFile.Data());
91   Int_t startTime = 0;
92   Int_t endTime   = 0;
93   for (Int_t run=0; run<nruns; run++){
94     Int_t irun=runArray[indexes[run]];
95     printf("Processing run %d ...\n",irun);
96     AliTPCcalibDB::Instance()->SetRun(irun);
97     dbutil->UpdateFromCalibDB();
98     //
99     AliDCSSensorArray *arrHV=calibDB->GetVoltageSensors(irun);
100     if (!arrHV) continue;
101     for  (Int_t isenHV=0; isenHV<arrHV->NumSensors(); ++isenHV){
102       AliDCSSensor *senHV=arrHV->GetSensorNum(isenHV);
103       if (!senHV) {
104         printf("Not interesting OCDB info\n");
105         continue;
106       }
107       startTime=senHV->GetStartTime();
108       endTime  =senHV->GetEndTime();
109       if (startTime>0&&endTime>0) break;
110     }
111     dbutil->FilterCE(120., 3., 4.,pcstream);
112     dbutil->FilterTracks(irun, 10.,pcstream);
113     AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(irun);          
114     //if (goofieArray) dbutil->FilterGoofie(goofieArray,2,4,pcstream);
115     // don't filter goofie for the moment
116     ProcessRun(irun, startTime,endTime);
117   }
118   delete pcstream;  
119 }
120
121
122 void ProcessRun(Int_t irun, Int_t startTime, Int_t endTime){
123   //
124   //
125   // 
126   AliSplineFit *fitVdrift=0x0;
127   Int_t startTimeGRP=0, stopTimeGRP=0;
128   if (calibDB->GetGRP(irun)){
129     startTimeGRP = AliTPCcalibDB::GetGRP(irun)->GetTimeStart();
130     stopTimeGRP  = AliTPCcalibDB::GetGRP(irun)->GetTimeEnd();
131   }
132   
133   AliTPCSensorTempArray * tempArray = AliTPCcalibDB::Instance()->GetTemperatureSensor(irun);
134   AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
135   AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(irun);
136   //
137   Int_t dtime = TMath::Max((endTime-startTime)/20,10*60);
138   //
139   //Goofie statistical data
140   //
141   TVectorD vecEntries, vecMean, vecMedian,vecRMS;
142   dbutil->ProcessGoofie(vecEntries ,vecMedian, vecMean, vecRMS);
143   //
144   //CE data processing - see ProcessCEdata function for description of the results
145   //
146   TVectorD fitResultsA, fitResultsC;
147   Int_t nmaskedCE;
148   Double_t chi2ACE=0,chi2CCE=0;
149   //     dbutil->ProcessCEdata("(sector<36)++gx++gy++lx++lx**2",fitResultsA,fitResultsC,nmaskedCE);
150   dbutil->ProcessCEdata("(sector<36)++gx++gy++(lx-134)++(sector<36)*(lx-134)++(ly/lx)^2",fitResultsA,fitResultsC,nmaskedCE,chi2ACE,chi2CCE);
151
152   TVectorD fitCEResultsA(7), fitCEResultsC(7);
153   Int_t    noutCE;
154   Double_t chi2CEA=0,chi2CEC=0;
155   AliTPCCalPad *time0 = dbutil->CreatePadTime0CE(fitCEResultsA, fitCEResultsC, noutCE, chi2CEA, chi2CEC);
156   delete time0;
157   //  
158   //
159   TVectorD vecTEntries, vecTMean, vecTRMS, vecTMedian, vecQEntries, vecQMean, vecQRMS, vecQMedian;
160   Float_t driftTimeA, driftTimeC;
161   dbutil->ProcessCEgraphs(vecTEntries, vecTMean, vecTRMS, vecTMedian,
162                           vecQEntries, vecQMean, vecQRMS, vecQMedian,
163                           driftTimeA, driftTimeC );
164   //
165   //
166   //
167   //drift velocity using tracks
168   //
169   //     fitVdrift=calibDB->GetVdriftSplineFit("ALISPLINEFIT_MEAN_VDRIFT_COSMICS_ALL",irun);
170   fitVdrift=calibDB->CreateVdriftSplineFit("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL",irun);
171   //noise data Processing - see ProcessNoiseData function for description of the results
172   TVectorD vNoiseMean, vNoiseMeanSenRegions, vNoiseRMS, vNoiseRMSSenRegions;
173   Int_t nonMaskedZero=0;
174   dbutil->ProcessNoiseData(vNoiseMean, vNoiseMeanSenRegions, vNoiseRMS, vNoiseRMSSenRegions, nonMaskedZero);
175   //
176   // comparisons
177   //
178   TVectorF pedestalDeviations;
179   TVectorF noiseDeviations;
180   TVectorF pulserQdeviations;
181   Float_t varQMean;
182   Int_t npadsOutOneTB;
183   Int_t npadsOffAdd;
184   dbutil->ProcessPedestalVariations(pedestalDeviations);
185   dbutil->ProcessNoiseVariations(noiseDeviations);
186   dbutil->ProcessPulserVariations(pulserQdeviations,varQMean,npadsOutOneTB,npadsOffAdd);
187   //
188   //L3 data 
189   //
190   Float_t bz=AliTPCcalibDB::GetBz(irun);
191   Char_t  l3pol=AliTPCcalibDB::GetL3Polarity(irun);
192   //
193   //calibration Pulser data processing
194   //
195   Int_t nOffChannels=0;
196   TVectorD vTimePulser;
197   nOffChannels=dbutil->GetNPulserOutliers();
198   dbutil->ProcessPulser(vTimePulser);
199   //
200   //ALTRO data
201   //
202   Int_t nMasked=0;
203   dbutil->ProcessALTROConfig(nMasked);
204   //
205   //Calib RAW data
206   //
207   Int_t nFailL1=-1;
208   if (calibDB->GetCalibRaw()) nFailL1=calibDB->GetCalibRaw()->GetNFailL1Phase();
209   //
210   //production information
211   //
212   Int_t nalien=0,nRawAlien=0,nlocal=0,nRawLocal=0;
213   //     GetProductionInfo(irun, nalien, nRawAlien, nlocal,nRawLocal);
214   //run type
215   TObjString runType(AliTPCcalibDB::GetRunType(irun).Data());
216   //
217   //
218   //
219   
220   for (Int_t itime=startTime; itime<endTime; itime+=dtime){
221     //
222     TTimeStamp tstamp(itime);
223     Float_t valuePressure  = calibDB->GetPressure(tstamp,irun,0);
224     Float_t valuePressure2 = calibDB->GetPressure(tstamp,irun,1);
225     Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative((UInt_t)itime,irun,0);
226     Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative((UInt_t)itime,irun,1);
227     //temperature fits
228     TLinearFitter * fitter = 0;
229     TVectorD vecTemp[10];
230     for (Int_t itype=0; itype<5; itype++)
231       for (Int_t iside=0; iside<2; iside++){
232         fitter= tempMap->GetLinearFitter(itype,iside,tstamp);     
233         if (!fitter) continue;
234         fitter->Eval();
235         fitter->GetParameters(vecTemp[itype+iside*5]);
236         delete fitter;
237       }
238     //
239     //measured skirt temperatures
240     //
241     TVectorD vecSkirtTempA(18);
242     TVectorD vecSkirtTempC(18);
243     Int_t nsenTemp=tempArray->NumSensors();
244     for (Int_t isenTemp=0;isenTemp<nsenTemp;++isenTemp){
245       AliTPCSensorTemp *senTemp=(AliTPCSensorTemp*)tempArray->GetSensorNum(isenTemp);
246       if (senTemp->GetType()!=3) continue;
247       if (TMath::Sqrt(senTemp->GetX()*senTemp->GetX()+senTemp->GetY()*senTemp->GetY())<100) continue; //only skirt, outer FC vessel
248       Double_t val=senTemp->GetValue(tstamp);
249       if (senTemp->GetSide()==0)
250         vecSkirtTempA[senTemp->GetSector()]=val;
251       else
252         vecSkirtTempC[senTemp->GetSector()]=val;
253     }
254     //
255     //goofie data
256     //
257     TVectorD vecGoofie;
258     if (goofieArray){
259       vecGoofie.ResizeTo(goofieArray->NumSensors());
260       for (Int_t isensor=0; isensor<goofieArray->NumSensors();isensor++){
261         AliDCSSensor *gsensor = goofieArray->GetSensor(isensor);
262         if (gsensor){
263           vecGoofie[isensor] = gsensor->GetValue(tstamp);
264         }
265       }
266     } else {
267       vecGoofie.ResizeTo(19);
268     }
269     //
270     TVectorD voltagesIROC(36);
271     TVectorD voltagesOROC(36);
272     for(Int_t j=1; j<36; j++) voltagesIROC[j-1] = AliTPCcalibDB::Instance()->GetChamberHighVoltage(irun, j,itime);
273     for(Int_t j=36; j<72; j++) voltagesOROC[j-36] = AliTPCcalibDB::Instance()->GetChamberHighVoltage(irun, j,itime);
274     Double_t voltIROC = TMath::Median(36, voltagesIROC.GetMatrixArray());
275     Double_t voltOROC = TMath::Median(36, voltagesOROC.GetMatrixArray());
276     //
277     Float_t  coverIA=AliTPCcalibDB::GetCoverVoltage(irun,0,itime);
278     Float_t  coverIC=AliTPCcalibDB::GetCoverVoltage(irun,18,itime);
279     Float_t  coverOA=AliTPCcalibDB::GetCoverVoltage(irun,36,itime);
280     Float_t  coverOC=AliTPCcalibDB::GetCoverVoltage(irun,54,itime);
281     Float_t  skirtA=AliTPCcalibDB::GetSkirtVoltage(irun,0,itime);
282     Float_t  skirtC=AliTPCcalibDB::GetSkirtVoltage(irun,18,itime);
283     Float_t  ggOffA=AliTPCcalibDB::GetGGoffsetVoltage(irun,0,itime);
284     Float_t  ggOffC=AliTPCcalibDB::GetGGoffsetVoltage(irun,18,itime);
285     //drift velocity
286     Float_t dvCorr=-5;
287     if (fitVdrift) dvCorr=fitVdrift->Eval(itime);
288     //
289     // Gain - Alexander Kalveit
290     //
291     Float_t  gainCosmic = 0;
292     Float_t  gainMIP = 0;
293     TObjArray * gainSplines = AliTPCcalibDB::Instance()->GetTimeGainSplinesRun(irun);
294     if (gainSplines) {
295       TGraphErrors * graphMIP = (TGraphErrors *) gainSplines->FindObject("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");
296       TGraphErrors * graphCosmic = (TGraphErrors *) gainSplines->FindObject("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL");
297       if (graphMIP) gainMIP = graphMIP->Eval(itime);
298       if (graphCosmic) gainCosmic = graphCosmic->Eval(itime);
299     }
300     
301     
302     //tempMap->GetLinearFitter(0,0,itime);
303     (*pcstream)<<"dcs"<<
304       "run="<<irun<<
305       "time="<<itime<<
306       "startTimeGRP="<<startTimeGRP<<
307       "stopTimeGRP="<<stopTimeGRP<<
308       //run type
309       "runType.="<<&runType<<
310       // voltage setting
311       "VIROC.="<<&voltagesIROC<<
312       "VOROC.="<<&voltagesOROC<<
313       "medianVIROC="<<voltIROC<<
314       "medianVOROC="<<voltOROC<<
315       "coverIA=" << coverIA <<
316       "coverIC=" << coverIC <<
317       "coverOA=" << coverOA <<
318       "coverOC=" << coverOC <<
319       "skirtA=" << skirtA <<
320       "skirtC=" << skirtC <<
321       "ggOffA=" << ggOffA <<
322       "ggOffC=" << ggOffC <<
323       //
324       "ptrel0="<<ptrelative0<<  // deltaTP/TP  - A side
325       "ptrel1="<<ptrelative1<<  // deltaTP/TPC - C side
326       "goofie.="<<&vecGoofie<<
327       "goofieE.="<<&vecEntries<<
328       "goofieMean.="<<&vecMean<<
329       "goofieMedian.="<<&vecMedian<<
330       "goofieRMS.="<<&vecRMS<<
331       //
332       "press="<<valuePressure<<
333       "press2="<<valuePressure2<<
334       "temp00.="<<&vecTemp[0]<<
335       "temp10.="<<&vecTemp[1]<<
336       "temp20.="<<&vecTemp[2]<<
337       "temp30.="<<&vecTemp[3]<<
338       "temp40.="<<&vecTemp[4]<<
339       "temp01.="<<&vecTemp[5]<<
340       "temp11.="<<&vecTemp[6]<<
341       "temp21.="<<&vecTemp[7]<<
342       "temp31.="<<&vecTemp[8]<<
343       "temp41.="<<&vecTemp[9]<<
344       "tempSkirtA.="<<&vecSkirtTempA<<
345       "tempSkirtC.="<<&vecSkirtTempC;
346     ProcessDrift(irun, itime);
347     ProcessDriftCE(irun,itime);
348     ProcessDriftAll(irun,itime);
349     ProcessKryptonTime(irun,itime);
350     (*pcstream)<<"dcs"<<        
351       //noise data
352       "meanNoise.="<<&vNoiseMean<<
353       "meanNoiseSen.="<<&vNoiseMeanSenRegions<<
354       "rmsNoise.="<<&vNoiseRMS<<
355       "rmsNoiseSen.="<<&vNoiseRMSSenRegions<<
356       "zeroNoise="<<nonMaskedZero<<
357       //pulser data
358       "timePulser.=" << &vTimePulser <<
359       "nOffPulser="<<nOffChannels<<
360       //altro data
361       "nMasked="<<nMasked<<
362       //ce data -Jens version
363       "CEfitA.="<<&fitResultsA<<
364       "CEfitC.="<<&fitResultsC<<
365       "nmaskedCE="<<nmaskedCE<<
366       "chi2ACE="<<chi2ACE<<
367       "chi2CCE="<<chi2CCE<<
368       //ce data new - MI version
369       "CEfitAMI.="<<&fitCEResultsA<<
370       "CEfitCMI.="<<&fitCEResultsC<<
371       "chi2CEA="<<chi2CEA<<
372       "chi2CEC="<<chi2CEC<<
373       //
374       //ce graph data
375       "CEgrTEntries.="<<&vecTEntries<<
376       "CEgrTMean.="<<&vecTMean<<
377       "CEgrTRMS.="<<&vecTRMS<<
378       "CEgrTMedian.="<<&vecTMedian<<
379       "CEgrQEntries.="<<&vecQEntries<<
380       "CEgrQMean.="<<&vecQMean<<
381       "CEgrQRMS.="<<&vecQRMS<<
382       "CEgrQMedian.="<<&vecQMedian<<
383       "CEgrDriftA="<<driftTimeA<<
384       "CEgrDriftC="<<driftTimeC<<
385       //calib raw data
386       "nFailL1="<<nFailL1<<
387       // b field
388       "Bz="<< bz <<
389       "L3polarity="<<l3pol<<
390       // production information
391       "nalien="<<nalien<<
392       "nRawAlien="<<nRawAlien<<
393       "nlocal="<<nlocal<<
394       "nRawLocal="<<nRawLocal<<
395       //comparisons with ref data
396       "pedestalDeviations.="<<&pedestalDeviations<<
397       "noiseDeviations.="<<&noiseDeviations<<
398       "pulserQdeviations.="<<&pulserQdeviations<<
399       //         "pulserVarQMean="<<varQMean<<
400       "pulserNpadsOutOneTB="<<npadsOutOneTB<<
401       "pulserNpadsOffAdd="<<npadsOffAdd<<
402       "driftCorrCosmAll="<<dvCorr<<
403       // time dependence of gain
404       "gainMIP="<<gainMIP<<
405       "gainCosmic="<<gainCosmic<<      
406       "\n";
407   }//end run loop
408   //     delete fitVdrift;
409   //     fitVdrift=0; 
410 }
411
412
413
414
415 void GetProductionInfo(Int_t run, Int_t &nalien, Int_t &nRawAlien, Int_t &nlocal, Int_t &nRawLocal){
416   //
417   // find number of ESDs in central and local production for run
418   //
419
420   nalien=0;
421   nRawAlien=0;
422   nlocal=0;
423   nRawLocal=0;
424   TString sNlines;
425   //find number of ESDs in alien
426   TString command="alien_find /alice/data/2009 ";
427   command += Form("%09d",run);
428   command += " | grep AliESDs.root | wc -l";
429   sNlines = gSystem->GetFromPipe(command.Data());
430   nalien=sNlines.Atoi();
431   //find number of raw files on alien
432   command="alien_find /alice/data/2009 ";
433   command += Form("%09d",run);
434   command += " | grep raw | grep -v tag | wc -l";
435   sNlines = gSystem->GetFromPipe(command.Data());
436   nRawAlien=sNlines.Atoi();
437   //find number of ESDs local
438   command="find /lustre/alice/alien/alice/data/2009 -name AliESDs.root | grep ";
439   command += Form("%09d",run);
440   command += " | wc -l";
441   sNlines = gSystem->GetFromPipe(command.Data());
442   nlocal=sNlines.Atoi();
443   //find number of local raw data files
444   command="find /lustre/alice/alien/alice/data/2009 -name \"*.root\" | grep ";
445   command += Form("%09d",run);
446   command += " | grep raw | grep -v tag | wc -l";
447   sNlines = gSystem->GetFromPipe(command.Data());
448   nRawLocal=sNlines.Atoi();
449 }
450
451
452
453 void ProcessDrift(Int_t run, Int_t timeStamp){
454   //
455   // dump drift calibration data to the tree
456   //
457   TObjArray *array =AliTPCcalibDB::Instance()->GetTimeVdriftSplineRun(run);
458   TGraphErrors *laserA[3]={0,0,0};
459   TGraphErrors *laserC[3]={0,0,0};
460   TGraphErrors *cosmicAll=0;
461   TGraphErrors *laserAE[3]={0,0,0};
462   TGraphErrors *laserCE[3]={0,0,0};
463   TGraphErrors *cosmicAllE=0;
464   static Double_t     vlaserA[3]={0,0,0};
465   static Double_t     vlaserC[3]={0,0,0};
466   static Double_t     vcosmicAll=0;
467   static Double_t     vdrift1=0;
468   vdrift1=AliTPCcalibDB::Instance()->GetVDriftCorrectionTime(timeStamp,run,0,1);
469
470   if (array){
471     laserA[0]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DELAY_LASER_ALL_A");
472     laserA[1]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_A");
473     laserA[2]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
474     laserC[0]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DELAY_LASER_ALL_C");
475     laserC[1]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_C");
476     laserC[2]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
477     cosmicAll =(TGraphErrors*)array->FindObject("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL");
478   }
479   if (laserA[0]) vlaserA[0]= AliTPCcalibDButil::EvalGraphConst(laserA[0],timeStamp);
480   if (laserA[1]) vlaserA[1]= AliTPCcalibDButil::EvalGraphConst(laserA[1],timeStamp);
481   if (laserA[2]) vlaserA[2]= AliTPCcalibDButil::EvalGraphConst(laserA[2],timeStamp);
482   if (laserC[0]) vlaserC[0]= AliTPCcalibDButil::EvalGraphConst(laserC[0],timeStamp);
483   if (laserC[1]) vlaserC[1]= AliTPCcalibDButil::EvalGraphConst(laserC[1],timeStamp);
484   if (laserC[2]) vlaserC[2]= AliTPCcalibDButil::EvalGraphConst(laserC[2],timeStamp);
485   if (cosmicAll) vcosmicAll= AliTPCcalibDButil::EvalGraphConst(cosmicAll,timeStamp); 
486   (*pcstream)<<"dcs"<<
487     "vlaserA0="<<vlaserA[0]<<
488     "vlaserA1="<<vlaserA[1]<<
489     "vlaserA2="<<vlaserA[2]<<
490     "vlaserC0="<<vlaserC[0]<<
491     "vlaserC1="<<vlaserC[1]<<
492     "vlaserC2="<<vlaserC[2]<<
493     "vcosmicAll="<<vcosmicAll<<
494     "vdrift1="<<vdrift1;
495
496   //
497   // define distance to measurement
498   //
499   static Double_t dlaserA=0; 
500   static Double_t dlaserC=0; 
501   static Double_t dcosmic=0; 
502   static Double_t slaserA=0; 
503   static Double_t slaserC=0; 
504   static Double_t scosmic=0; 
505   static Double_t  vclaserA[3]={0,0,0};
506   static Double_t  vclaserC[3]={0,0,0};
507   static Double_t  vccosmicAll=0;
508   for (Int_t i=0;i<3;i++){
509     if (laserA[i]) AliTPCcalibDButil::GetNearest(laserA[i],timeStamp,dlaserA,vclaserA[i]);
510     if (laserC[i]) AliTPCcalibDButil::GetNearest(laserC[i],timeStamp,dlaserC,vclaserC[i]);
511   }  
512   if (cosmicAll) AliTPCcalibDButil::GetNearest(cosmicAll,timeStamp,dcosmic,vccosmicAll);
513   (*pcstream)<<"dcs"<<
514     "vclaserA0="<<vclaserA[0]<<
515     "vclaserA1="<<vclaserA[1]<<
516     "vclaserA2="<<vclaserA[2]<<
517     "vclaserC0="<<vclaserC[0]<<
518     "vclaserC1="<<vclaserC[1]<<
519     "vclaserC2="<<vclaserC[2]<<
520     "vccosmicAll="<<vccosmicAll<<
521     "dlaserA="<<dlaserA<<
522     "dlaserC="<<dlaserC<<
523     "dcosmic="<<dcosmic<<
524     "slaserA="<<slaserA<<
525     "slaserC="<<slaserC<<
526     "scosmic="<<scosmic;
527 }
528
529 void ProcessDriftCE(Int_t run,Int_t timeStamp){
530   //
531   // dump drift calibration data CE
532   //
533   TObjArray *arrT=AliTPCcalibDB::Instance()->GetCErocTtime();
534   AliTPCParam *param=AliTPCcalibDB::Instance()->GetParameters();
535   static TVectorD tdriftCE(74);
536   static TVectorD tndriftCE(74);
537   static TVectorD vdriftCE(74);
538   static TVectorD tcdriftCE(74);
539   static TVectorD tddriftCE(74);
540   static Double_t ltime0A;
541   static Double_t ltime0C;
542   //
543   //
544   //
545   ltime0A  = dbutil->GetLaserTime0(run,timeStamp,36000,0);
546   ltime0C  = dbutil->GetLaserTime0(run,timeStamp,36000,1);
547   //
548   for (Int_t i=0; i<arrT->GetEntries();i++){
549     tdriftCE[i]=0;
550     vdriftCE[i]=0;
551     TGraph *graph = (TGraph*)arrT->At(i);
552     if (!graph) continue;
553     tdriftCE[i]=AliTPCcalibDButil::EvalGraphConst(graph,timeStamp);
554     Double_t deltaT,gry;
555     AliTPCcalibDButil::GetNearest(graph,timeStamp,deltaT,gry);
556     tndriftCE[i]=graph->GetN();
557     tcdriftCE[i]=gry;          
558     tddriftCE[i]=deltaT;               
559     if (i%36<18){
560       vdriftCE[i] =(param->GetZLength(i)/(tdriftCE[i]*param->GetTSample()*(1.-ltime0A)-param->GetL1Delay()))/param->GetDriftV();
561     }else{
562       vdriftCE[i] =(param->GetZLength(i)/(tdriftCE[i]*param->GetTSample()*(1.-ltime0A)-param->GetL1Delay()))/param->GetDriftV();
563     }
564   }
565
566   (*pcstream)<<"dcs"<<  
567     "tdriftCE.="<<&tdriftCE<<
568     "vdriftCE.="<<&vdriftCE<<
569     "tndriftCE.="<<&tndriftCE<<
570     "tcdriftCE.="<<&tcdriftCE<<
571     "tddriftCE.="<<&tddriftCE<<
572     "ltime0A="<<ltime0A<<
573     "ltime0C="<<ltime0C;
574    }
575
576
577 void ProcessDriftAll(Int_t run,Int_t timeStamp){
578   //
579   // dump drift calibration data  all calibrations form DB util
580   // test of utils
581   static Double_t vdriftCEA=0, vdriftCEC=0, vdriftCEM=0;
582   static Double_t vdriftLTA=0, vdriftLTC=0, vdriftLTM=0;
583   static Double_t dce=0, dla=0,dlc=0,dlm=0;
584   static Double_t ltime0A;
585   static Double_t ltime0C;
586   static Double_t     vdrift1=0;
587   vdrift1=AliTPCcalibDB::Instance()->GetVDriftCorrectionTime(timeStamp,run,0,1);
588   
589   //
590   vdriftCEA= dbutil->GetVDriftTPCCE(dla,run,timeStamp,36000,0);
591   vdriftCEC= dbutil->GetVDriftTPCCE(dlc,run,timeStamp,36000,1);
592   vdriftCEM= dbutil->GetVDriftTPCCE(dlm,run,timeStamp,36000,2);
593   //
594   vdriftLTA= dbutil->GetVDriftTPCLaserTracks(dla,run,timeStamp,36000,0);
595   vdriftLTC= dbutil->GetVDriftTPCLaserTracks(dlc,run,timeStamp,36000,1);
596   vdriftLTM= dbutil->GetVDriftTPCLaserTracks(dlm,run,timeStamp,36000,2);
597   //
598   ltime0A  = dbutil->GetLaserTime0(run,timeStamp,36000,0);
599   ltime0C  = dbutil->GetLaserTime0(run,timeStamp,36000,1);
600
601   (*pcstream)<<"dcs"<<  
602     //
603     "vdriftCEA="<<vdriftCEA<<
604     "vdriftCEC="<<vdriftCEC<<
605     "vdriftCEM="<<vdriftCEM<<
606     "dce="<<dce<<
607     "vdriftLTA="<<vdriftLTA<<
608     "vdriftLTC="<<vdriftLTC<<
609     "vdriftLTM="<<vdriftLTM<<
610     "dla="<<dla<<
611     "dlc="<<dlc<<
612     "dlm="<<dlm<<
613     "vdrift1="<<vdrift1;
614
615 }
616
617
618
619 void ProcessKryptonTime(Int_t run, Int_t timeStamp){
620   //
621   // Dumping  krypton calibration results
622   //
623   static TObjArray * krArray=0;
624   if (!krArray) {
625     AliCDBEntry* entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeGainKrypton", run);
626     if (entry){
627       krArray = (TObjArray*)entry->GetObject();
628     }
629   }
630   static TVectorD krMean(74);
631   static TVectorD krErr(74);
632   static TVectorD krDist(74);
633   TGraphErrors *gr=0;
634   Double_t deltaT=0,gry=0;
635   if (krArray){
636     for (Int_t isec=0; isec<72; isec++){
637       krMean[isec]=0;
638       krDist[isec]=0;
639       krErr[isec]=0;
640       gr=(TGraphErrors*)krArray->At(isec);
641       if (gr) {
642         krMean[isec]=gr->Eval(timeStamp);
643         AliTPCcalibDButil::GetNearest(gr, timeStamp,deltaT,gry);
644         krDist[isec]=deltaT;
645       }     
646       if (72+isec<krArray->GetEntries()) {
647         gr=(TGraphErrors*)krArray->At(72+isec);
648         if (gr) krErr[isec]=gr->Eval(timeStamp);
649       }
650     }
651     krMean[72]= TMath::Median(36,krMean.GetMatrixArray());
652     krMean[73]= TMath::Median(36,&(krMean.GetMatrixArray()[36]));
653     krErr[72]= TMath::Median(36,krErr.GetMatrixArray());
654     krErr[73]= TMath::Median(36,&(krErr.GetMatrixArray()[36]));
655   }
656   (*pcstream)<<"dcs"<<
657     "krMean.="<<&krMean<<
658     "krErr.="<<&krErr<<
659     "krDist.="<<&krDist;
660 }
661
662
663
664