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