]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/TPCcalib/AliTPCcalibSummary.cxx
TPC module
[u/mrichter/AliRoot.git] / TPC / TPCcalib / 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 #include "AliTPCCalibGlobalMisalignment.h"
62 #include "AliTPCExBTwist.h"
63 #include "AliTPCComposedCorrection.h"
64 //
65 //
66 //
67 AliTPCcalibSummary::AliTPCcalibSummary():
68   TNamed(),
69   fCalibDB(0),
70   fDButil(0),
71   fPcstream(0)
72 {
73   //
74   // default constructor
75   // OCDB have to be setupe before - not part of class
76   // usualy ConfigOCDB.C macro used
77   // 
78   fPcstream = new TTreeSRedirector("dcsTime.root");
79   fCalibDB = AliTPCcalibDB::Instance();
80   fDButil= new AliTPCcalibDButil;
81 }
82
83 AliTPCcalibSummary::~AliTPCcalibSummary(){
84   //
85   // destructor  - close streamer
86   //
87   delete fPcstream;  
88 }
89
90 void AliTPCcalibSummary::Process(const char * runList, Int_t first, Int_t last){
91   //
92   // runList - listOfRuns to process
93   // first   - first run to process
94   // last    - last  to process
95   // 
96   //
97   // make list of runs
98   //
99
100   ifstream inputFile;
101   inputFile.open("run.list");
102   Int_t irun=0;
103   TArrayI runArray(100000);
104   Int_t indexes[100000];
105   Int_t nruns=0;
106   printf("Runs to process:\n");
107   if (!inputFile.is_open()) { 
108     printf("Problem to open file %s\n",runList);
109   }
110   while( inputFile.good() ) {
111     inputFile >> irun;
112     printf("Run \t%d\n",irun);
113     if (irun<first) continue;  // process only subset of list
114     if (last>0 && irun>last) continue;  // process only subset of list
115     runArray[nruns]=irun;
116     nruns++;
117   }
118
119
120   TMath::Sort(nruns, runArray.fArray, indexes,kFALSE);
121   Int_t startTime = 0;
122   Int_t endTime   = 0;
123   for (Int_t run=0; run<nruns; run++){
124     irun=runArray[indexes[run]];
125     printf("Processing run %d ...\n",irun);
126     fCalibDB->SetRun(irun);
127     fDButil->UpdateFromCalibDB();
128     fDButil->SetReferenceRun(irun);
129     fDButil->UpdateRefDataFromOCDB();
130     fCalibDB->CreateGUITree("calPads.root");
131     fDButil->CreateGUIRefTree("calPadsRef.root");
132     //
133     AliDCSSensorArray *arrHV=fCalibDB->GetVoltageSensors(irun);
134     if (!arrHV) continue;
135     for  (Int_t isenHV=0; isenHV<arrHV->NumSensors(); ++isenHV){
136       AliDCSSensor *senHV=arrHV->GetSensorNum(isenHV);
137       if (!senHV) {
138         printf("Not interesting OCDB info\n");
139         continue;
140       }
141       startTime=senHV->GetStartTime();
142       endTime  =senHV->GetEndTime();
143       if (startTime>0&&endTime>0) break;
144     }
145     AliDCSSensorArray* goofieArray = fCalibDB->GetGoofieSensors(irun);  
146     if (goofieArray) fDButil->FilterGoofie(goofieArray,0.5,4.,4,10,fPcstream);
147     // don't filter goofie for the moment
148     ProcessRun(irun, startTime,endTime);
149   }
150 }
151
152
153 void AliTPCcalibSummary::ProcessRun(Int_t irun, Int_t startTime, Int_t endTime){
154   //
155   // Process run irun
156   // 
157   fCalibDB->SetRun(irun);
158   fDButil->UpdateFromCalibDB();
159   fDButil->SetReferenceRun(irun);
160   fDButil->UpdateRefDataFromOCDB();
161   fCalibDB->CreateGUITree("calPads.root");
162   fDButil->CreateGUIRefTree("calPadsRef.root");
163
164   //
165   AliSplineFit *fitVdrift=0x0;
166   Int_t startTimeGRP=0, stopTimeGRP=0;
167   if (fCalibDB->GetGRP(irun)){
168     startTimeGRP = AliTPCcalibDB::GetGRP(irun)->GetTimeStart();
169     stopTimeGRP  = AliTPCcalibDB::GetGRP(irun)->GetTimeEnd();
170   }
171   if (startTime==0){
172     startTime=startTimeGRP;
173     endTime=stopTimeGRP;
174   }
175   AliTPCSensorTempArray * tempArray = fCalibDB->GetTemperatureSensor(irun);
176   AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
177   AliDCSSensorArray* goofieArray = fCalibDB->GetGoofieSensors(irun);
178   //
179   Int_t dtime = TMath::Max((endTime-startTime)/20,10);
180   //
181   //Goofie statistical data
182   //
183   TVectorD vecEntries, vecMean, vecMedian,vecRMS;
184   fDButil->ProcessGoofie(vecEntries ,vecMedian, vecMean, vecRMS);
185   //
186   //CE data processing - see ProcessCEdata function for description of the results
187   //
188   TVectorD fitResultsA, fitResultsC;
189   Int_t nmaskedCE;
190   Double_t chi2ACE=0,chi2CCE=0;
191   //     fDButil->ProcessCEdata("(sector<36)++gx++gy++lx++lx**2",fitResultsA,fitResultsC,nmaskedCE);
192   fDButil->ProcessCEdata("(sector<36)++gx++gy++(lx-134)++(sector<36)*(lx-134)++(ly/lx)^2",fitResultsA,fitResultsC,nmaskedCE,chi2ACE,chi2CCE);
193
194   TVectorD fitCEResultsA(7), fitCEResultsC(7);
195   Int_t    noutCE;
196   Double_t chi2CEA=0,chi2CEC=0;
197   AliTPCCalPad *time0 = fDButil->CreatePadTime0CE(fitCEResultsA, fitCEResultsC, noutCE, chi2CEA, chi2CEC);
198   delete time0;
199   //  
200   //
201   TVectorD vecTEntries, vecTMean, vecTRMS, vecTMedian, vecQEntries, vecQMean, vecQRMS, vecQMedian;
202   Float_t driftTimeA, driftTimeC;
203   fDButil->ProcessCEgraphs(vecTEntries, vecTMean, vecTRMS, vecTMedian,
204                           vecQEntries, vecQMean, vecQRMS, vecQMedian,
205                           driftTimeA, driftTimeC );
206   //
207   //
208   //
209   //drift velocity using tracks
210   //
211   //     fitVdrift=fCalibDB->GetVdriftSplineFit("ALISPLINEFIT_MEAN_VDRIFT_COSMICS_ALL",irun);
212   fitVdrift=fCalibDB->CreateVdriftSplineFit("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL",irun);
213   //noise data Processing - see ProcessNoiseData function for description of the results
214   TVectorD vNoiseMean, vNoiseMeanSenRegions, vNoiseRMS, vNoiseRMSSenRegions;
215   Int_t nonMaskedZero=0, nNaN=0;
216   fDButil->ProcessNoiseData(vNoiseMean, vNoiseMeanSenRegions, vNoiseRMS, vNoiseRMSSenRegions, nonMaskedZero, nNaN);
217   //
218   // comparisons
219   //
220   TVectorF pedestalDeviations;
221   TVectorF noiseDeviations;
222   TVectorF pulserQdeviations;
223   Float_t varQMean;
224   Int_t npadsOutOneTB;
225   Int_t npadsOffAdd;
226   fDButil->ProcessPedestalVariations(pedestalDeviations);
227   fDButil->ProcessNoiseVariations(noiseDeviations);
228   fDButil->ProcessPulserVariations(pulserQdeviations,varQMean,npadsOutOneTB,npadsOffAdd);
229   //
230   //L3 data 
231   //
232   Float_t bz=AliTPCcalibDB::GetBz(irun);
233   Char_t  l3pol=AliTPCcalibDB::GetL3Polarity(irun);
234   //
235   //QA data processing
236   //
237   TVectorD vQaOcc;
238   TVectorD vQaQtot;
239   TVectorD vQaQmax;
240   fDButil->ProcessQAData(vQaOcc, vQaQtot, vQaQmax);
241   //
242   //calibration Pulser data processing
243   //
244   Int_t nOffChannels=0;
245   TVectorD vTimePulser;
246   nOffChannels=fDButil->GetNPulserOutliers();
247   fDButil->ProcessPulser(vTimePulser);
248   //
249   //ALTRO data
250   //
251   Int_t nMasked=0;
252   fDButil->ProcessALTROConfig(nMasked);
253   //
254   //Calib RAW data
255   //
256   Int_t nFailL1=-1;
257   if (fCalibDB->GetCalibRaw()) nFailL1=fCalibDB->GetCalibRaw()->GetNFailL1Phase();
258   //
259   //production information
260   //
261   Int_t nalien=0,nRawAlien=0,nlocal=0,nRawLocal=0;
262   //run type
263   TObjString runType(AliTPCcalibDB::GetRunType(irun).Data());
264   //
265   //
266   //
267   
268   for (Int_t itime=startTime; itime<endTime; itime+=dtime){
269     //
270     TTimeStamp tstamp(itime);
271     Float_t valuePressure  = fCalibDB->GetPressure(tstamp,irun,0);
272     Float_t valuePressure2 = fCalibDB->GetPressure(tstamp,irun,1);
273     Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative((UInt_t)itime,irun,0);
274     Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative((UInt_t)itime,irun,1);
275     //temperature fits
276     TLinearFitter * fitter = 0;
277     TVectorD vecTemp[10];
278     for (Int_t itype=0; itype<5; itype++)
279       for (Int_t iside=0; iside<2; iside++){
280         fitter= tempMap->GetLinearFitter(itype,iside,tstamp);
281         if (!fitter) continue;
282         fitter->Eval();
283         fitter->GetParameters(vecTemp[itype+iside*5]);
284         delete fitter;
285       }
286     //
287     //measured skirt temperatures
288     //
289     TVectorD vecSkirtTempA(18);
290     TVectorD vecSkirtTempC(18);
291     Int_t nsenTemp=tempArray->NumSensors();
292     for (Int_t isenTemp=0;isenTemp<nsenTemp;++isenTemp){
293       AliTPCSensorTemp *senTemp=(AliTPCSensorTemp*)tempArray->GetSensorNum(isenTemp);
294       if (senTemp->GetType()!=3) continue;
295       if (TMath::Sqrt(senTemp->GetX()*senTemp->GetX()+senTemp->GetY()*senTemp->GetY())<100) continue; //only skirt, outer FC vessel
296       Double_t val=senTemp->GetValue(tstamp);
297       if (senTemp->GetSide()==0)
298         vecSkirtTempA[senTemp->GetSector()]=val;
299       else
300         vecSkirtTempC[senTemp->GetSector()]=val;
301     }
302     //
303     //goofie data
304     //
305     TVectorD vecGoofie; 
306     if (goofieArray){
307       vecGoofie.ResizeTo(goofieArray->NumSensors());
308       for (Int_t isensor=0; isensor<goofieArray->NumSensors();isensor++){
309         AliDCSSensor *gsensor = goofieArray->GetSensor(isensor);
310         if (gsensor){
311           vecGoofie[isensor] = gsensor->GetValue(tstamp);
312         }
313       }
314     } else {
315       vecGoofie.ResizeTo(19);
316     }
317     //
318     static TVectorF voltagesIROC(36);
319     static TVectorF voltagesIROCMedian(36);
320     static TVectorF voltagesIROCNominal(36);
321     static TVectorF voltagesIROCCurrentNominal(36);
322     static TVectorF voltagesIROCStatus(36);
323     static TVectorF voltagesIROCGoodFraction(36);
324     //
325     static TVectorF voltagesOROC(36);
326     static TVectorF voltagesOROCMedian(36);
327     static TVectorF voltagesOROCNominal(36);
328     static TVectorF voltagesOROCCurrentNominal(36);
329     static TVectorF voltagesOROCStatus(36);
330     static TVectorF voltagesOROCGoodFraction(36);
331     
332     for(Int_t j=0; j<36; j++){
333       voltagesIROC[j]               = fCalibDB->GetChamberHighVoltage(irun, j,itime);
334       voltagesIROCMedian[j]         = fCalibDB->GetChamberHighVoltageMedian(j);
335       voltagesIROCNominal[j]        = fCalibDB->GetParameters()->GetNominalVoltage(j);
336       voltagesIROCCurrentNominal[j] = fCalibDB->GetChamberCurrentNominalHighVoltage(j);
337       voltagesIROCStatus[j]         = fCalibDB->GetChamberHVStatus(j);
338       voltagesIROCGoodFraction[j]   = fCalibDB->GetChamberGoodHighVoltageFraction(j);
339     }
340     
341     for(Int_t j=36; j<72; j++) {
342       voltagesOROC[j-36]               = fCalibDB->GetChamberHighVoltage(irun, j,itime);
343       voltagesOROCMedian[j-36]         = fCalibDB->GetChamberHighVoltageMedian(j);
344       voltagesOROCNominal[j-36]        = fCalibDB->GetParameters()->GetNominalVoltage(j);
345       voltagesOROCCurrentNominal[j-36] = fCalibDB->GetChamberCurrentNominalHighVoltage(j);
346       voltagesOROCStatus[j-36]         = fCalibDB->GetChamberHVStatus(j);
347       voltagesOROCGoodFraction[j-36]   = fCalibDB->GetChamberGoodHighVoltageFraction(j);
348     }
349     
350     Double_t voltIROC = TMath::Median(36, voltagesIROC.GetMatrixArray());
351     Double_t voltOROC = TMath::Median(36, voltagesOROC.GetMatrixArray());
352     //
353     Float_t  coverIA=AliTPCcalibDB::GetCoverVoltage(irun,0,itime);
354     Float_t  coverIC=AliTPCcalibDB::GetCoverVoltage(irun,18,itime);
355     Float_t  coverOA=AliTPCcalibDB::GetCoverVoltage(irun,36,itime);
356     Float_t  coverOC=AliTPCcalibDB::GetCoverVoltage(irun,54,itime);
357     Float_t  skirtA=AliTPCcalibDB::GetSkirtVoltage(irun,0,itime);
358     Float_t  skirtC=AliTPCcalibDB::GetSkirtVoltage(irun,18,itime);
359     Float_t  ggOffA=AliTPCcalibDB::GetGGoffsetVoltage(irun,0,itime);
360     Float_t  ggOffC=AliTPCcalibDB::GetGGoffsetVoltage(irun,18,itime);
361     //drift velocity
362     Float_t dvCorr=-5;
363     if (fitVdrift) dvCorr=fitVdrift->Eval(itime);
364     //data taking active
365     Bool_t dataTakingActive=fCalibDB->IsDataTakingActive((time_t)itime);
366     
367     //tempMap->GetLinearFitter(0,0,itime);
368     // parameter description is as follows: // <name>;<unit>;<group>;<array type/array info>
369     //      <name> should also include the unit; e.g. IROC anode voltage (V))
370     // where group is e.g.
371     //       o HV:              for High Voltage information
372     //       o Environment:     for environmental information (temperature; pressure)
373     //       o Pulser:          Calibartion pulser information
374     //       o CE:              Central electrode inforamtion
375     //       o Noise Pedestals: Noise and Pedestal information
376     //       o ALTRO:           ALTRO configuration information
377     //
378     // <array type> describes the information stored in the elements of a TVectorT<> variable
379     // this can e.g. be
380     //       o Sector:   The array has 36 entries; one per sector (A00-A17 and C00 to C17)
381     //       o Sector-A: The array has 18 entries; one per sector on the A-Side (A00-A17)
382     //       o Sector-C: The array has 18 entries; one per sector on the C-Side (C00-C17)
383     //       o ROC:      The array has 72 entries; one per ROC (IA00-IA17, IC00-IC17, OA00-OA17, OC00-OC17)
384     //       o Parameters from a fit: in this case the description per parameter should be given
385     (*fPcstream)<<"dcs"<<
386       "run="<<irun<< // Run number
387       "time="<<itime<< // Time stamp of calibration entry
388       "startTimeGRP="<<startTimeGRP<< // Start time of run from GRP
389       "stopTimeGRP="<<stopTimeGRP<< // Stop time of run from GRP
390       "dataTakingActive="<<dataTakingActive<< // If data taking is active
391       //run type
392       "runType.="<<&runType<< // Run Type; e.g. PHYSICS; LASER; COSMIC; PEDESTAL; PULSER
393       // voltage setting
394       "VIROC.="               << &voltagesIROC               << // IROC anode voltage [calib interval] (V);HV;Sector
395       "VIROCMedian.="         << &voltagesIROCMedian         << // IROC anode voltage [Median of run] (V);HV;Sector
396       "VIROCNominal.="        << &voltagesIROCNominal        << // IROC anode voltage [global nominal] (V);HV;Sector
397       "VIROCCurrentNominal.=" << &voltagesIROCCurrentNominal << // IROC anode voltage [current nominal] (V);HV;Sector
398       "VIROCGoodHVFraction.=" << &voltagesIROCGoodFraction   << // IROC anode voltage [fraction of good settings];-;HV;Sector
399       "VIROCStatus.="         << &voltagesIROCStatus         << // IROC HV status;-;HV;Sector
400       //
401       "VOROC.="               << &voltagesOROC               << // OROC anode voltage [calib interval] (V);HV;Sector
402       "VOROCMedian.="         << &voltagesOROCMedian         << // OROC anode voltage [Median of run] (V);HV;Sector
403       "VOROCNominal.="        << &voltagesOROCNominal        << // OROC anode voltage [global nominal] (V);HV;Sector
404       "VOROCCurrentNominal.=" << &voltagesOROCCurrentNominal << // OROC anode voltage [current nominal] (V);HV;Sector
405       "VOROCGoodHVFraction.=" << &voltagesOROCGoodFraction   << // OROC anode voltage [fraction of good settings];-;HV;Sector
406       "VOROCStatus.="         << &voltagesOROCStatus         << // OROC HV status;-;HV;Sector
407       //
408       "medianVIROC="          << voltIROC                    << // IROC anode voltage [median of all IROCs] (V);HV
409       "medianVOROC="          << voltOROC                    << // OROC anode voltage [median of all OROCs] (V);HV
410       "coverIA="              << coverIA                     << // Cover voltage IROC A-Side (V);HV
411       "coverIC="              << coverIC                     << // Cover voltage IROC C-Side (V);HV
412       "coverOA="              << coverOA                     << // Cover voltage OROC A-Side (V);HV
413       "coverOC="              << coverOC                     << // Cover voltage OROC C-Side (V);HV
414       "skirtA="               << skirtA                      << // Skirt voltage IROC A-Side (V);HV
415       "skirtC="               << skirtC                      << // Skirt voltage IROC C-Side (V);HV
416       "ggOffA="               << ggOffA                      << // Gating grid offset voltage A-Side (V);HV
417       "ggOffC="               << ggOffC                      << // Gating grid offset voltage C-Side (V);HV
418       //
419       "ptrel0="<<ptrelative0<<  // deltaTP/TP  - A side
420       "ptrel1="<<ptrelative1<<  // deltaTP/TPC - C side
421       "goofie.="<<&vecGoofie<<
422       "goofieE.="<<&vecEntries<<
423       "goofieMean.="<<&vecMean<<
424       "goofieMedian.="<<&vecMedian<<
425       "goofieRMS.="<<&vecRMS<<
426       //
427       "press="<<valuePressure<<
428       "press2="<<valuePressure2<<
429       "temp00.="<<&vecTemp[0]<< // T-Fit ROC A-Side;Environment;Mean Temp (#circC);dT/dgx (K/cm);dT/dgy (K/cm)
430       "temp10.="<<&vecTemp[1]<< // T-Fit OFC A-Side;Environment;Mean Temp (#circC);dT/dz (K/cm);dT/d#phi (K/rad)
431       "temp20.="<<&vecTemp[2]<< // T-Fit IFC+TS A-Side;Environment;Mean Temp (#circC);dT/dz (K/cm);dT/d#phi (K/rad)
432       "temp30.="<<&vecTemp[3]<< // T-Fit Skirt A-Side;Environment;Mean Temp (#circC);dT/dgx (K/cm);dT/dgy (K/cm)
433       "temp40.="<<&vecTemp[4]<< // T-Fit IFC A-Side;Environment;Mean Temp (#circC);dT/dz (K/cm);dT/d#phi (K/rad)
434       "temp01.="<<&vecTemp[5]<< // T-Fit ROC C-Side;Environment;Mean Temp (#circC);dT/dgx (K/cm);dT/dgy (K/cm)
435       "temp11.="<<&vecTemp[6]<< // T-Fit OFC C-Side;Environment;Mean Temp (#circC);dT/dz (K/cm);dT/d#phi (K/rad)
436       "temp21.="<<&vecTemp[7]<< // T-Fit IFC+TS C-Side;Environment;Mean Temp (#circC);dT/dz (K/cm);dT/d#phi (K/rad)
437       "temp31.="<<&vecTemp[8]<< // T-Fit Skirt C-Side;Environment;Mean Temp (#circC);dT/dgx (K/cm);dT/dgy (K/cm)
438       "temp41.="<<&vecTemp[9]<< // T-Fit IFC C-Side;Environment;Mean Temp (#circC);dT/dz (K/cm);dT/d#phi (K/rad)
439       //
440       "tempSkirtA.="<<&vecSkirtTempA<< // T Skirt A-Side;Environment;Sector-A
441       "tempSkirtC.="<<&vecSkirtTempC;  // T Skirt C-Side;Environment;Sector-C
442
443     ProcessDrift(irun, itime);
444     ProcessDriftCE(irun,itime);
445     ProcessDriftAll(irun,itime);
446     //    ProcessKryptonTime(irun,itime);
447     ProcessCTP(irun,itime);
448     ProcessAlign(irun,itime);
449     ProcessGain(irun,itime);
450     //ProcessDriftCERef();
451     //ProcessPulserRef();
452     //ProcessCurrent(irun,itime);
453
454
455     (*fPcstream)<<"dcs"<<       
456       //noise data
457       "meanNoise.="    << &vNoiseMean           << // Mean Noise;Noise Pedestals;All Pads;IROCs;OROCs small pads;OROCs large pads
458       "meanNoiseSen.=" << &vNoiseMeanSenRegions << // Mean Noise in sensitive regions;Noise Pedestals;All Pads;IROCs;OROCs small pads;OROCs large pads
459       "rmsNoise.="     << &vNoiseRMS            << // RMS Noise;Noise Pedestals;All Pads;IROCs;OROCs small pads;OROCs large pads
460       "rmsNoiseSen.="  << &vNoiseRMSSenRegions  << // RMS Noise in sensitive regions;Noise Pedestals;All Pads;IROCs;OROCs small pads;OROCs large pads
461       "zeroNoise="     << nonMaskedZero         << // Pads with zero noise;Noise Pedestals
462       "nNaN="          << nNaN                  << // Pads with NaN noise;Noise Pedestals
463       //QA data
464       "occQA.="  << &vQaOcc  <<
465       "qQA.="    << &vQaQtot <<
466       "qmaxQA.=" << &vQaQmax <<
467       //pulser data
468       "timePulser.=" << &vTimePulser <<
469       "nOffPulser="<<nOffChannels<<
470       //altro data
471       "nMasked="<< nMasked << // Number of masked pads;ALTRO
472       //
473       //ce data -Jens version
474       //
475       "CEfitA.="        << &fitResultsA  << // CE-Fit A-Side;CE;Offset (timebins);IROC/OROC Offset (timebins);dt/dgx (timebins/cm);dt/dgy (timebins/cm);dt/dlx common (timebins/cm);dt/dlx IROCs (timebins/cm)
476       "CEfitC.="        << &fitResultsC  << // CE Fit C-Side;CE;Offset (timebins);IROC/OROC Offset (timebins);dt/dgx (timebins/cm);dt/dgy (timebins/cm);dt/dlx common (timebins/cm);dt/dlx IROCs (timebins/cm)
477       "nmaskedCE="      << nmaskedCE     << // CE Number of outliers;CE
478       "chi2ACE="        << chi2ACE       << // CE-Fit Chi^{2} A-Side;CE
479       "chi2CCE="        << chi2CCE       << // CE-Fit Chi^{2} C-Side;CE
480       //
481       //ce data new - MI version
482       //
483       "CEfitAMI.="      << &fitCEResultsA<< // CE-Fit A-Side [MI];CE;Offset (timebins);IROC/OROC Offset (timebins);dt/dgx (timebins/cm);dt/dgy (timebins/cm);dt/dlx common (timebins/cm);dt/dlx IROCs (timebins/cm)
484       "CEfitCMI.="      << &fitCEResultsC<< // CE-Fit C-Side [MI];CE;Offset (timebins);IROC/OROC Offset (timebins);dt/dgx (timebins/cm);dt/dgy (timebins/cm);dt/dlx common (timebins/cm);dt/dlx IROCs (timebins/cm)
485       "chi2CEA="        << chi2CEA       << // CE-Fit Chi^{2} A-Side [MI];CE
486       "chi2CEC="        << chi2CEC       << // CE-Fit Chi^{2} C-Side [MI];CE
487       //
488       //ce graph data
489       //
490       "CEgrTEntries.="   << &vecTEntries << // CE-graph drift time - entries;CE;ROC
491       "CEgrTMean.="      << &vecTMean    << // CE-graph mean drift time;CE;ROC
492       "CEgrTRMS.="       << &vecTRMS     << // CE-graph RMS of drift time;CE;ROC
493       "CEgrTMedian.="    << &vecTMedian  << // CE-graph median drift time;CE;ROC
494       "CEgrQEntries.="   << &vecQEntries << // CE-graph charge - entries;CE;ROC
495       "CEgrQMean.="      << &vecQMean    << // CE-graph mean charge;CE;ROC
496       "CEgrQRMS.="       << &vecQRMS     << // CE-graph RMS charge;CE;ROC
497       "CEgrQMedian.="    << &vecQMedian  << // CE-graph median charge;CE;ROC
498       "CEgrDriftA="      << driftTimeA   << // CE median drift time A-Side;CE
499       "CEgrDriftC="      << driftTimeC   << // CE median drift time C-Side;CE
500       //
501       //calib raw data
502       //
503       "nFailL1="         << nFailL1      << // RCU synchonisation failures;ALTRO
504       // b field
505       "Bz="              << bz           << // Magnetic Field (T);Environment
506       "L3polarity="      << l3pol        << // L3 polarity;Environment
507       // production information
508       "nalien="          << nalien       << // obsolete
509       "nRawAlien="       << nRawAlien    << // obsolete
510       "nlocal="          << nlocal       << // obsolete
511       "nRawLocal="       <<nRawLocal     << // obsolete
512       //
513       // comparisons with ref data
514       //
515       "pedestalDeviations.=" << &pedestalDeviations << // Pedestal variation to ref (fraction);Noise Pedestals;>#pm 0.5 ADC;>#pm 1 ADC;>#pm 1.5 ADC;>#pm 2.0 ADC
516       "noiseDeviations.="    << &noiseDeviations    << // Noise var to ref (fraction);Noise Pedestals;>5%;>10%;>15%;>20%
517       "pulserQdeviations.="  << &pulserQdeviations  << // Pulser-Q var to ref (fraction);Pulser;>0.5%;>1%;>5%;>10%
518       //         "pulserVarQMean="<<varQMean<<
519       "pulserNpadsOutOneTB=" << npadsOutOneTB       << // Number of pads with Pulser time var >#pm 1 tb to ROC mean;Pulser
520       "pulserNpadsOffAdd="   << npadsOffAdd         << // Number of pads without signal but signal in ref;Pulser
521       "driftCorrCosmAll="    << dvCorr              <<
522       "\n";
523   }//end run loop
524 }
525
526
527
528
529
530
531 void AliTPCcalibSummary::ProcessDrift(Int_t run, Int_t timeStamp){
532   //
533   // dump drift calibration data to the tree
534   //
535   TObjArray *array =fCalibDB->GetTimeVdriftSplineRun(run);
536   TGraphErrors *laserA[3]={0,0,0};
537   TGraphErrors *laserC[3]={0,0,0};
538   TGraphErrors *cosmicAll=0;
539   static Double_t     vlaserA[3]={0,0,0};
540   static Double_t     vlaserC[3]={0,0,0};
541   static Double_t     vcosmicAll=0;
542   static Double_t     vdrift1=0;                                // TODO: repeated below, obsolete?
543   vdrift1=fCalibDB->GetVDriftCorrectionTime(timeStamp,run,0,1); // TODO: repeated below, obsolete?
544
545   if (array){
546     laserA[0]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DELAY_LASER_ALL_A");
547     laserA[1]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_A");
548     laserA[2]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
549     laserC[0]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DELAY_LASER_ALL_C");
550     laserC[1]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_C");
551     laserC[2]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
552     cosmicAll =(TGraphErrors*)array->FindObject("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL");
553   }
554
555   //
556   // TODO: the information stored in vlaserXX, vcosmicAll and vclaserXX,vccosmicAll 
557   //       seem to be redundant information do we need to keep vlaserXX vcosmicAll
558   //
559   if (laserA[0]) vlaserA[0]= AliTPCcalibDButil::EvalGraphConst(laserA[0],timeStamp);
560   if (laserA[1]) vlaserA[1]= AliTPCcalibDButil::EvalGraphConst(laserA[1],timeStamp);
561   if (laserA[2]) vlaserA[2]= AliTPCcalibDButil::EvalGraphConst(laserA[2],timeStamp);
562   if (laserC[0]) vlaserC[0]= AliTPCcalibDButil::EvalGraphConst(laserC[0],timeStamp);
563   if (laserC[1]) vlaserC[1]= AliTPCcalibDButil::EvalGraphConst(laserC[1],timeStamp);
564   if (laserC[2]) vlaserC[2]= AliTPCcalibDButil::EvalGraphConst(laserC[2],timeStamp);
565   if (cosmicAll) vcosmicAll= AliTPCcalibDButil::EvalGraphConst(cosmicAll,timeStamp);
566   (*fPcstream)<<"dcs"<<
567     "vlaserA0="   << vlaserA[0] << // Laser offset A-Side;Drift            //TODO: Obsolete
568     "vlaserA1="   << vlaserA[1] << // Laser drift correction A-Side;Drift  //TODO: Obsolete
569     "vlaserA2="   << vlaserA[2] << // Laser gy correction A-Side;Drift     //TODO: Obsolete
570     "vlaserC0="   << vlaserC[0] << // Laser offset C-Side;Drift            //TODO: Obsolete
571     "vlaserC1="   << vlaserC[1] << // Laser drift correction C-Side;Drift  //TODO: Obsolete
572     "vlaserC2="   << vlaserC[2] << // Laser gy correction C-Side;Drift     //TODO: Obsolete
573     "vcosmicAll=" << vcosmicAll << // Cosmic drift corrrection;Drift       //TODO: Obsolete
574     //
575     "vdrift1="    << vdrift1;      // Combined drift correction ;Drift      // TODO: repeated below, obsolete?
576
577   //
578   // define distance to measurement
579   //
580   static Double_t dlaserA=0; 
581   static Double_t dlaserC=0; 
582   static Double_t dcosmic=0; 
583   static Double_t slaserA=0; //TODO: Obsolete?
584   static Double_t slaserC=0; //TODO: Obsolete?
585   static Double_t scosmic=0; //TODO: Obsolete?
586   static Double_t  vclaserA[3]={0,0,0};
587   static Double_t  vclaserC[3]={0,0,0};
588   static Double_t  vccosmicAll=0;
589   for (Int_t i=0;i<3;i++){
590     if (laserA[i]) AliTPCcalibDButil::GetNearest(laserA[i],timeStamp,dlaserA,vclaserA[i]);
591     if (laserC[i]) AliTPCcalibDButil::GetNearest(laserC[i],timeStamp,dlaserC,vclaserC[i]);
592   }  
593   if (cosmicAll) AliTPCcalibDButil::GetNearest(cosmicAll,timeStamp,dcosmic,vccosmicAll);
594   (*fPcstream)<<"dcs"<<
595     "vclaserA0="   << vclaserA[0]<< // Laser offset A-Side;Drift
596     "vclaserA1="   << vclaserA[1]<< // Laser drift correction A-Side;Drift
597     "vclaserA2="   << vclaserA[2]<< // Laser gy correction A-Side;Drift
598     "vclaserC0="   << vclaserC[0]<< // Laser offset A-Side;Drift
599     "vclaserC1="   << vclaserC[1]<< // Laser drift correction A-Side;Drift
600     "vclaserC2="   << vclaserC[2]<< // Laser gy correction A-Side;Drift
601     "vccosmicAll=" << vccosmicAll<< // Cosmic drift corrrection;Drift
602     "dlaserA="     << dlaserA    << // Distance to laser measurement A-Side;Drift
603     "dlaserC="     << dlaserC    << // Distance to laser measurement C-Side;Drift
604     "dcosmic="     << dcosmic    << // Distance to cosmics measurement A-Side;Drift
605     "slaserA="     << slaserA    << //TODO: Obsolete?
606     "slaserC="     << slaserC    << //TODO: Obsolete?
607     "scosmic="     << scosmic;      //TODO: Obsolete?
608
609   static TGeoMatrix * matrixAlign=0;
610   static Double_t twistX=0;
611   static Double_t twistY=0;
612   if (matrixAlign==0){
613     AliTPCComposedCorrection * corr =  (AliTPCComposedCorrection *)array->FindObject("FitCorrectionTime");
614     if (!corr) {
615       matrixAlign=new TGeoHMatrix;
616       
617     }
618     if (corr){
619        AliTPCCalibGlobalMisalignment *align = (AliTPCCalibGlobalMisalignment*)corr->GetCorrections()->FindObject("FitAlignTime");
620        AliTPCExBTwist *twist  = (AliTPCExBTwist*)corr->GetCorrections()->FindObject("FitExBTwistTime");
621        if (twist){
622          twistX=twist->GetXTwist();
623          twistY=twist->GetYTwist();
624          //delete twist;
625        }
626        if (align && align->GetAlignGlobal()){
627          matrixAlign =  (TGeoMatrix*) (align->GetAlignGlobal()->Clone());        
628          //delete align;
629        }
630     }    
631   }
632   (*fPcstream)<<"dcs"<<
633     "alignTime.="<<matrixAlign<<
634     "twistX="<<twistX<<
635     "twistY="<<twistY;
636 }
637
638 void AliTPCcalibSummary::ProcessDriftCE(Int_t run,Int_t timeStamp){
639   //
640   // dump drift calibration data CE
641   //
642   TObjArray *arrT=fCalibDB->GetCErocTtime();
643   AliTPCParam *param=fCalibDB->GetParameters();
644   static TVectorD tdriftCE(74);
645   static TVectorD tndriftCE(74);
646   static TVectorD vdriftCE(74);
647   static TVectorD tcdriftCE(74);
648   static TVectorD tddriftCE(74);
649   static Double_t ltime0A=0.;
650   static Double_t ltime0C=0.;
651   //
652   //
653   //
654   ltime0A  = fDButil->GetLaserTime0(run,timeStamp,36000,0);
655   ltime0C  = fDButil->GetLaserTime0(run,timeStamp,36000,1);
656   //
657   // TODO: are the values tdriftCE and tcdriftCE redundant and need both be kept?
658   for (Int_t i=0; i<arrT->GetEntries();i++){
659     tdriftCE[i]=0;
660     vdriftCE[i]=0;
661     TGraph *graph = (TGraph*)arrT->At(i);
662     if (!graph) continue;
663     tdriftCE[i]=AliTPCcalibDButil::EvalGraphConst(graph,timeStamp);
664     Double_t deltaT,gry;
665     AliTPCcalibDButil::GetNearest(graph,timeStamp,deltaT,gry);
666     tndriftCE[i]=graph->GetN();
667     tcdriftCE[i]=gry;
668     tddriftCE[i]=deltaT;
669     // TODO: use tdriftCE or tndriftCE
670     if (i%36<18){
671       vdriftCE[i] =(param->GetZLength(i)/(tdriftCE[i]*param->GetTSample()*(1.-ltime0A)-param->GetL1Delay()))/param->GetDriftV();
672     }else{
673       vdriftCE[i] =(param->GetZLength(i)/(tdriftCE[i]*param->GetTSample()*(1.-ltime0A)-param->GetL1Delay()))/param->GetDriftV();
674     }
675   }
676   // export values
677   (*fPcstream)<<"dcs"<<  
678     "tdriftCE.="  << &tdriftCE  << // CE arrival time;CE;Sector                  // TODO: obsolete, redundant?
679     "vdriftCE.="  << &vdriftCE  << // CE derived drift velocity;CE;Sector
680     "tndriftCE.=" << &tndriftCE << // CE number of points;CE;Sector
681     "tcdriftCE.=" << &tcdriftCE << // CE arrival time - nearest point;CE;Sector
682     "tddriftCE.=" << &tddriftCE << // CE distance to closest measuement;CE;Sector
683     "ltime0A="    << ltime0A    << // CE laser offset A-Side;CE
684     "ltime0C="    << ltime0C     ; // CE laser offset C-Side;CE
685    }
686
687
688 void AliTPCcalibSummary::ProcessDriftAll(Int_t run,Int_t timeStamp){
689   //
690   // dump drift calibration data  all calibrations form DB util
691   // test of utils
692   static Double_t vdriftCEA=0, vdriftCEC=0, vdriftCEM=0;
693   static Double_t vdriftLTA=0, vdriftLTC=0, vdriftLTM=0;
694   static Double_t vdriftLTAon=0, vdriftLTCon=0, vdriftLTMon=0;
695   static Double_t vdriftITS=0;
696   static Double_t vdriftP=0;
697   static Double_t dcea=0, dcec=0, dcem=0,  dla=0,dlc=0,dlm=0, dlaon=0,dlcon=0,dlmon=0, dp=0;
698   static Double_t dits=0;
699   static Double_t ltime0A=0.;
700   static Double_t ltime0C=0.;
701   static Double_t ctime0=0.;
702   static Double_t vdrift1=0;
703   vdrift1= fCalibDB->GetVDriftCorrectionTime(timeStamp,run,0,1);
704   vdriftP = fDButil->GetVDriftTPC(dp, run, timeStamp, 86400, 3600,0);
705   ctime0 = AliTPCcalibDButil::GetTriggerOffsetTPC(run,timeStamp, 36000, 3600,0);
706   //
707   vdriftCEA= fDButil->GetVDriftTPCCE(dcea,run,timeStamp,36000,0);
708   vdriftCEC= fDButil->GetVDriftTPCCE(dcec,run,timeStamp,36000,1);
709   vdriftCEM= fDButil->GetVDriftTPCCE(dcem,run,timeStamp,36000,2);
710   //
711   vdriftLTA= fDButil->GetVDriftTPCLaserTracks(dla,run,timeStamp,36000,0);
712   vdriftLTC= fDButil->GetVDriftTPCLaserTracks(dlc,run,timeStamp,36000,1);
713   vdriftLTM= fDButil->GetVDriftTPCLaserTracks(dlm,run,timeStamp,36000,2);
714   //
715   vdriftLTAon= fDButil->GetVDriftTPCLaserTracksOnline(dlaon,run,timeStamp,36000,0);
716   vdriftLTCon= fDButil->GetVDriftTPCLaserTracksOnline(dlcon,run,timeStamp,36000,1);
717   vdriftLTMon= fDButil->GetVDriftTPCLaserTracksOnline(dlmon,run,timeStamp,36000,2);
718   //
719   vdriftITS= fDButil->GetVDriftTPCITS(dits, run,timeStamp);
720   //
721   ltime0A  = fDButil->GetLaserTime0(run,timeStamp,36000,0); // TODO: not used, needed?
722   ltime0C  = fDButil->GetLaserTime0(run,timeStamp,36000,1); // TODO: not used, needed?
723
724   (*fPcstream)<<"dcs"<<  
725     //
726     "vdriftCEA="   << vdriftCEA   << // CE drift correction A-Side;Drift
727     "vdriftCEC="   << vdriftCEC   << // CE drift correction C-Side;Drift
728     "vdriftCEM="   << vdriftCEM   << // CE drift correction Mean;Drift
729     "dcea="        << dcea        << // CE distance to closest measurement A-Side;Drift
730     "dcec="        << dcec        << // CE distance to closest measurement C-Side;Drift
731     "dcem="        << dcem        << // CE distance to closest measurement Mean;Drift
732     "vdriftLTA="   << vdriftLTA   << // Offline Laser track vdrift correction A-Side;Drift
733     "vdriftLTC="   << vdriftLTC   << // Offline Laser track vdrift correction C-Side;Drift
734     "vdriftLTM="   << vdriftLTM   << // Offline Laser track vdrift correction Mean;Drift
735     "dla="         << dla         << // Offline Laser track distance to closest measurement A-Side;Drift
736     "dlc="         << dlc         << // Offline Laser track distance to closest measurement C-Side;Drift
737     "dlm="         << dlm         << // Offline Laser track distance to closest measurement Mean;Drift
738     "vdriftLTAon=" << vdriftLTAon << // Online Laser track vdrift correction A-Side;Drift
739     "vdriftLTCon=" << vdriftLTCon << // Online Laser track vdrift correction C-Side;Drift
740     "vdriftLTMon=" << vdriftLTMon << // Online Laser track vdrift correction Mean;Drift
741     "dlaOn="       << dlaon       << // Online Laser track distance to closest measurement A-Side;Drift
742     "dlcOn="       << dlcon       << // Online Laser track distance to closest measurement C-Side;Drift
743     "dlmOn="       << dlmon       << // Online Laser track distance to closest measurement Mean;Drift
744     //
745     //
746     "vdriftITS="   << vdriftITS   << // TPC-ITS vdrift correction;Drift
747     "dits="        << dits        << // TPC-ITS vdrift correction distance to closest measurement;Drift
748     "ctime0="      << ctime0      << // Trigger offset correction;Drift
749     "vdriftP="     << vdriftP     << // Cosmics vdrift correction;Drift
750     "dp="          << dp          << // Cosmics vdrift correction distance to closest measurement;Drift
751     "vdrift1="     << vdrift1      ; // combined drift velocity;Drift
752
753 }
754
755
756
757 void AliTPCcalibSummary::ProcessKryptonTime(Int_t run, Int_t timeStamp){
758   //
759   // Dumping  krypton calibration results
760   //
761   static TObjArray * krArray=0;
762   if (!krArray) {
763     AliCDBEntry* entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeGainKrypton", run);
764     if (entry){
765       krArray = (TObjArray*)entry->GetObject();
766     }
767   }
768   static TVectorD krMean(74);
769   static TVectorD krErr(74);
770   static TVectorD krDist(74);
771   TGraphErrors *gr=0;
772   Double_t deltaT=0,gry=0;
773   if (krArray){
774     for (Int_t isec=0; isec<72; isec++){
775       krMean[isec]=0;
776       krDist[isec]=0;
777       krErr[isec]=0;
778       gr=(TGraphErrors*)krArray->At(isec);
779       if (gr) {
780         krMean[isec]=AliTPCcalibDButil::EvalGraphConst(gr,timeStamp);
781         AliTPCcalibDButil::GetNearest(gr, timeStamp,deltaT,gry);
782         krDist[isec]=deltaT;
783       }     
784       if (72+isec<krArray->GetEntries()) {
785         gr=(TGraphErrors*)krArray->At(72+isec);
786         if (gr) krErr[isec]=AliTPCcalibDButil::EvalGraphConst(gr,timeStamp);
787       }
788     }
789     krMean[72]= TMath::Median(36,krMean.GetMatrixArray());
790     krMean[73]= TMath::Median(36,&(krMean.GetMatrixArray()[36]));
791     krErr[72]= TMath::Median(36,krErr.GetMatrixArray());
792     krErr[73]= TMath::Median(36,&(krErr.GetMatrixArray()[36]));
793   }
794   (*fPcstream)<<"dcs"<<
795     "krMean.="<<&krMean<<
796     "krErr.="<<&krErr<<
797     "krDist.="<<&krDist;
798 }
799
800 void AliTPCcalibSummary::ProcessCTP(Int_t irun, Int_t timeStamp){
801   //
802   // 
803   //
804   static TClonesArray *pcarray = new TClonesArray("AliCTPInputTimeParams",1);
805   static AliCTPTimeParams *ctpParams =0;
806   ctpParams = fCalibDB->GetCTPTimeParams(); //
807   const TObjArray        *arr = ctpParams->GetInputTimeParams();
808   pcarray->ExpandCreateFast(TMath::Max(arr->GetEntries(),1));
809   for (Int_t i=0; i<arr->GetEntries(); i++){
810     new ((*pcarray)[i]) AliCTPInputTimeParams(*((AliCTPInputTimeParams*)arr->At(i)));
811   }
812   (*fPcstream)<<"ctp"<<
813     "run="<<irun<<
814     "time="<<timeStamp<<
815     "ctpP.="<<ctpParams<<
816     "ctpArr="<<pcarray<<
817     "\n";
818   (*fPcstream)<<"dcs"<<
819     "ctpP.="<<ctpParams<<
820     "ctpArr="<<pcarray;
821 }
822
823 void AliTPCcalibSummary::ProcessGain(Int_t irun, Int_t timeStamp){
824   //
825   // Dump gain related information to the tree
826   //
827   static Float_t  gainCosmic = 0;
828   static Float_t  gainMIP = 0;
829   static Float_t  attachMIP = 0;
830   static Double_t dMIP=0; 
831   Double_t dummy=0;
832   static TVectorD vGainGraphIROC(36);
833   static TVectorD vGainGraphOROCmed(36);
834   static TVectorD vGainGraphOROClong(36);
835   static TVectorD vGainGraphIROCErr(36);
836   static TVectorD vGainGraphOROCmedErr(36);
837   static TVectorD vGainGraphOROClongErr(36);
838   //
839   static TVectorD vGainQMaxGraphRegion(3);
840   static TVectorD vGainQTotGraphRegion(3);
841   //
842   static TGraphErrors ggrPadEqualMax(36);
843   static TGraphErrors ggrPadEqualTot(36);
844   //
845   static TGraphErrors ggrDipAngleMaxShort;
846   static TGraphErrors ggrDipAngleMaxMedium;
847   static TGraphErrors ggrDipAngleMaxLong;
848   static TGraphErrors ggrDipAngleMaxAbsolute;
849   //
850   static TGraphErrors ggrDipAngleTotShort;
851   static TGraphErrors ggrDipAngleTotMedium;
852   static TGraphErrors ggrDipAngleTotLong;
853   static TGraphErrors ggrDipAngleTotAbsolute;
854   //
855   static TVectorD vFitDipAngleParMaxShort(3);
856   static TVectorD vFitDipAngleParMaxMedium(3);
857   static TVectorD vFitDipAngleParMaxLong(3);
858   static TVectorD vFitDipAngleParMaxAbsolute(3);
859   //
860   static TVectorD vFitDipAngleParTotShort(3);
861   static TVectorD vFitDipAngleParTotMedium(3);
862   static TVectorD vFitDipAngleParTotLong(3);
863   static TVectorD vFitDipAngleParTotAbsolute(3);
864
865   
866   vGainGraphIROC.Zero();
867   vGainGraphOROCmed.Zero();
868   vGainGraphOROClong.Zero();
869   vGainGraphIROCErr.Zero();
870   vGainGraphOROCmedErr.Zero();
871   vGainGraphOROClongErr.Zero();
872   vGainQMaxGraphRegion.Zero();
873   vGainQTotGraphRegion.Zero();
874   TGraphErrors grDummy;
875   TObjArray * gainSplines = fCalibDB->GetTimeGainSplinesRun(irun);
876   if (gainSplines) {
877     TGraphErrors * graphMIP = (TGraphErrors *) gainSplines->FindObject("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");
878     TGraphErrors * graphCosmic = (TGraphErrors *) gainSplines->FindObject("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL");
879     TGraphErrors * graphAttach = (TGraphErrors *) gainSplines->FindObject("TGRAPHERRORS_MEAN_ATTACHMENT_BEAM_ALL");
880     //
881     TGraphErrors * grPadEqualQMax = (TGraphErrors * ) gainSplines->FindObject("TGRAPHERRORS_MEANQMAX_PADREGIONGAIN_BEAM_ALL");
882     TGraphErrors * grPadEqualQTot = (TGraphErrors * ) gainSplines->FindObject("TGRAPHERRORS_MEANQTOT_PADREGIONGAIN_BEAM_ALL");
883     //
884     TGraphErrors * graphGainIROC       = (TGraphErrors *) gainSplines->FindObject("TGRAPHERRORS_MEAN_CHAMBERGAIN_SHORT_BEAM_ALL");
885     TGraphErrors * graphGainOROCMedium = (TGraphErrors *) gainSplines->FindObject("TGRAPHERRORS_MEAN_CHAMBERGAIN_MEDIUM_BEAM_ALL");
886     TGraphErrors * graphGainOROCLong   = (TGraphErrors *) gainSplines->FindObject("TGRAPHERRORS_MEAN_CHAMBERGAIN_LONG_BEAM_ALL");
887     //
888     //
889     TF1*  funDipAngleMax[4]={0x0,0x0,0x0,0x0};
890     TF1*  funDipAngleTot[4]={0x0,0x0,0x0,0x0};
891     TGraphErrors*  grDipAngleMax[4]={0x0,0x0,0x0,0x0};
892     TGraphErrors*  grDipAngleTot[4]={0x0,0x0,0x0,0x0};
893     const char* names[4]={"SHORT","MEDIUM","LONG","ABSOLUTE"};
894     for (Int_t iPadRegion=0; iPadRegion<4; ++iPadRegion) {
895       funDipAngleMax[iPadRegion]=(TF1*) gainSplines->FindObject(Form("TF1_QMAX_DIPANGLE_%s_BEAM_ALL",names[iPadRegion]));
896       funDipAngleTot[iPadRegion]=(TF1*) gainSplines->FindObject(Form("TF1_QTOT_DIPANGLE_%s_BEAM_ALL",names[iPadRegion]));
897       grDipAngleMax[iPadRegion]= (TGraphErrors*) gainSplines->FindObject(Form("TGRAPHERRORS_QMAX_DIPANGLE_%s_BEAM_ALL",names[iPadRegion]));
898       grDipAngleTot[iPadRegion]= (TGraphErrors*) gainSplines->FindObject(Form("TGRAPHERRORS_QTOT_DIPANGLE_%s_BEAM_ALL",names[iPadRegion]));
899     }
900     //
901     for(Int_t iPar=0; iPar < 3; iPar++) {
902       if (funDipAngleMax[0]) vFitDipAngleParMaxShort(iPar)    = funDipAngleMax[0]->GetParameter(iPar);
903       if (funDipAngleMax[1]) vFitDipAngleParMaxMedium(iPar)   = funDipAngleMax[1]->GetParameter(iPar);
904       if (funDipAngleMax[2]) vFitDipAngleParMaxLong(iPar)     = funDipAngleMax[2]->GetParameter(iPar);
905       if (funDipAngleMax[3]) vFitDipAngleParMaxAbsolute(iPar) = funDipAngleMax[3]->GetParameter(iPar);
906       //
907       if (funDipAngleTot[0]) vFitDipAngleParTotShort(iPar)    = funDipAngleTot[0]->GetParameter(iPar);
908       if (funDipAngleTot[1]) vFitDipAngleParTotMedium(iPar)   = funDipAngleTot[1]->GetParameter(iPar);
909       if (funDipAngleTot[2]) vFitDipAngleParTotLong(iPar)     = funDipAngleTot[2]->GetParameter(iPar);
910       if (funDipAngleTot[3]) vFitDipAngleParTotAbsolute(iPar) = funDipAngleTot[3]->GetParameter(iPar);
911     }
912     //
913     if (grDipAngleMax[0]) ggrDipAngleMaxShort    = * grDipAngleMax[0];
914     if (grDipAngleMax[1]) ggrDipAngleMaxMedium   = * grDipAngleMax[1];
915     if (grDipAngleMax[2]) ggrDipAngleMaxLong     = * grDipAngleMax[2];
916     if (grDipAngleMax[3]) ggrDipAngleMaxAbsolute = * grDipAngleMax[3];
917     //
918     if (grDipAngleTot[0]) ggrDipAngleTotShort    = * grDipAngleTot[0];
919     if (grDipAngleTot[1]) ggrDipAngleTotMedium   = * grDipAngleTot[1];
920     if (grDipAngleTot[2]) ggrDipAngleTotLong     = * grDipAngleTot[2];
921     if (grDipAngleTot[3]) ggrDipAngleTotAbsolute = * grDipAngleTot[3];
922     //
923     //
924     TGraphErrors *grPadEqualMax = (TGraphErrors * ) gainSplines->FindObject("TGRAPHERRORS_MEANQMAX_PADREGIONGAIN_BEAM_ALL");
925     TGraphErrors *grPadEqualTot = (TGraphErrors * ) gainSplines->FindObject("TGRAPHERRORS_MEANQTOT_PADREGIONGAIN_BEAM_ALL");
926     if (grPadEqualMax) ggrPadEqualMax = *grPadEqualMax;
927     if (grPadEqualTot) ggrPadEqualTot = *grPadEqualTot;
928
929
930     if (graphGainIROC && graphGainOROCMedium && graphGainOROCLong) {
931       Double_t x=0,y=0;
932       for (Int_t i=0; i<36; ++i){
933         graphGainIROC->GetPoint(i,x,y);
934         vGainGraphIROC(i)=y;
935         graphGainOROCMedium->GetPoint(i,x,y);
936         vGainGraphOROCmed(i)=y;
937         graphGainOROCLong->GetPoint(i,x,y);
938         vGainGraphOROClong(i)=y;
939         //errors
940         vGainGraphIROCErr(i)     = graphGainIROC->GetEY()[i];
941         vGainGraphOROCmedErr(i)  = graphGainOROCMedium->GetEY()[i];
942         vGainGraphOROClongErr(i) = graphGainOROCLong->GetEY()[i];
943       }
944        for (Int_t i=0; i<3; ++i){
945          vGainQMaxGraphRegion[i]=grPadEqualQMax->GetY()[i];
946          vGainQTotGraphRegion[i]=grPadEqualQTot->GetY()[i];
947        }
948     }
949     
950     if (graphMIP) gainMIP = AliTPCcalibDButil::EvalGraphConst(graphMIP,timeStamp);
951     if (graphCosmic) gainCosmic = AliTPCcalibDButil::EvalGraphConst(graphCosmic,timeStamp);
952     if (graphAttach) attachMIP = AliTPCcalibDButil::EvalGraphConst(graphAttach,timeStamp);
953     if (graphMIP)  AliTPCcalibDButil::GetNearest(graphMIP, timeStamp, dMIP,dummy);
954   }
955     
956   // time dependence of gain 
957   (*fPcstream)<<"dcs"<<
958     "grPadEqualMax.="             << &ggrPadEqualMax             <<
959     "grPadEqualTot.="             << &ggrPadEqualTot             <<
960     "rocGainIROC.="               << &vGainGraphIROC             <<
961     "rocGainOROCMedium.="         << &vGainGraphOROCmed          <<
962     "rocGainOROCLong.="           << &vGainGraphOROClong         <<
963     "rocGainErrIROC.="            << &vGainGraphIROCErr          <<
964     "rocGainErrOROCMedium.="      << &vGainGraphOROCmedErr       <<
965     "rocGainErrOROCLong.="        << &vGainGraphOROClongErr      <<
966     "vGainQMaxGraphRegion.="      << &vGainQMaxGraphRegion       <<
967     "vGainQTotGraphRegion.="      << &vGainQTotGraphRegion       <<
968     //
969     "vFitDipAngleParMaxShort.="   << &vFitDipAngleParMaxShort    <<
970     "vFitDipAngleParMaxMedium.="  << &vFitDipAngleParMaxMedium   <<
971     "vFitDipAngleParMaxLong.="    << &vFitDipAngleParMaxLong     <<
972     "vFitDipAngleParMaxAbsolute.="<< &vFitDipAngleParMaxAbsolute <<
973     //
974     "vFitDipAngleParTotShort.="   << &vFitDipAngleParTotShort    <<
975     "vFitDipAngleParTotMedium.="  << &vFitDipAngleParTotMedium   <<
976     "vFitDipAngleParTotLong.="    << &vFitDipAngleParTotLong     <<
977     "vFitDipAngleParTotAbsolute.="<< &vFitDipAngleParTotAbsolute <<
978     //
979     "grDipAngleMaxShort.="        << &ggrDipAngleMaxShort        <<
980     "grDipAngleMaxMedium.="       << &ggrDipAngleMaxMedium       <<
981     "grDipAngleMaxLong.="         << &ggrDipAngleMaxLong         <<
982     "grDipAngleMaxAbsolute.="     << &ggrDipAngleMaxAbsolute     <<
983     //
984     "grDipAngleTotShort.="        << &ggrDipAngleTotShort        <<
985     "grDipAngleTotMedium.="       << &ggrDipAngleTotMedium       <<
986     "grDipAngleTotLong.="         << &ggrDipAngleTotLong         <<
987     "grDipAngleTotAbsolute.="     << &ggrDipAngleTotAbsolute     <<
988     //
989     "gainMIP="                    << gainMIP                     <<
990     "attachMIP="                  << attachMIP                   <<
991     "dMIP="                       << dMIP                        <<
992     "gainCosmic="                 << gainCosmic                   ;
993 }
994
995
996 void AliTPCcalibSummary::ProcessAlign(Int_t run, Int_t timeStamp){
997   //
998   // Proccess alignment 
999   //
1000   TString   grnames[12]={"ALIGN_ITS", "ALIGN_ITSP", "ALIGN_ITSM", "ALIGN_ITSB",
1001                          "ALIGN_TRD", "ALIGN_TRDP", "ALIGN_TRDM","ALIGN_TRDB",
1002                          "ALIGN_TOF", "ALIGN_TOFP", "ALIGN_TOFM","ALIGN_TOFB"};
1003   TString   grpar[9]={"DELTAPSI", "DELTATHETA", "DELTAPHI",
1004                       "DELTAX", "DELTAY", "DELTAZ",
1005                       "DRIFTVD", "T0", "VDGY"};
1006   static Double_t values[12*9];
1007   static Double_t errs[12*9];
1008
1009   TObjArray *array =fCalibDB->GetTimeVdriftSplineRun(run);
1010   TGraphErrors *gr=0;
1011   for (Int_t idet=0; idet<12; idet++){
1012     for (Int_t ipar=0; ipar<9; ipar++){
1013       TString grName=grnames[idet];
1014       grName+="_TPC_";
1015       grName+=grpar[ipar];
1016       if (array){
1017         gr = (TGraphErrors*)array->FindObject(grName.Data());
1018       }
1019       values[9*idet+ipar]=0;
1020       errs[9*idet+ipar]=0;
1021       if (gr) values[9*idet+ipar]=AliTPCcalibDButil::EvalGraphConst(gr,timeStamp);
1022       (*fPcstream)<<"dcs"<<
1023         Form("%s=",grName.Data())<<values[9*idet+ipar];
1024       (*fPcstream)<<"align"<<
1025         Form("%s=",grName.Data())<<values[9*idet+ipar];
1026     }
1027   }
1028   (*fPcstream)<<"align"<<
1029     "time="<<timeStamp<<
1030     "run="<<run<<
1031     "\n";
1032 }
1033
1034
1035 void AliTPCcalibSummary::ProcessDriftCERef(){
1036   //
1037   // Get fit of residuals if CE in respect with reference
1038   // data
1039   //
1040   static TVectorD  sec(72);
1041   static TVectorD vec0(72);
1042   static TVectorD vecLy(72);
1043   static TVectorD vecLx(72);
1044   static TVectorD vecChi2(72);
1045   static TVectorD vecN(72);
1046   //
1047   static TVectorD vecA0(72);
1048   static TVectorD vecALy(72);
1049   static TVectorD vecALx(72);
1050   static TVectorD vecAChi2(72);
1051   //
1052   static TVectorD vecASide(4);
1053   static TVectorD vecCSide(4);
1054   static Bool_t isCalc=kFALSE;
1055   
1056   TFile f("calPads.root");
1057   TFile fref("calPadsRef.root");
1058   TTree * tree = (TTree*)f.Get("calPads");
1059   TTree * treeRef = (TTree*)fref.Get("calPads");
1060   tree->AddFriend(treeRef,"R");
1061   tree->SetAlias("inCE","((CEQmean.fElements>35)&&abs(CETmean.fElements)<1.5&&abs(CETrms.fElements/1.2-1)<0.2)");  // outlyerTrms
1062   tree->SetAlias("inCER","((R.CEQmean.fElements>35)&&abs(R.CETmean.fElements)<1.5&&abs(R.CETrms.fElements/1.2-1)<0.2)");  // outlyerTrms
1063   //
1064   if (!isCalc){
1065     // make fits only once
1066     TStatToolkit toolkit;
1067     Double_t chi2=0;
1068     Int_t    npoints=0;
1069     TVectorD param;
1070     TMatrixD covar;
1071     tree->SetAlias("dt","CETmean.fElements-R.CETmean.fElements");
1072     TCut cutAll ="inCE&&inCER&&abs(CETmean.fElements-R.CETmean.fElements)<0.5";
1073     TString  fstringG="";              // global part
1074     fstringG+="ly.fElements++";
1075     fstringG+="(lx.fElements-134.)++";
1076     for (Int_t isec=0; isec<72; isec++){
1077       TStatToolkit::FitPlane(tree,"0.264*dt", fstringG.Data(),Form("sector==%d",isec)+cutAll, chi2,npoints,param,covar,-1,0, 10000000, kFALSE);
1078       if (npoints<3) continue;
1079       printf("Sector=%d\n",isec);
1080       vec0[isec]=param[0];
1081       vecLy[isec]=param[1];
1082       vecLx[isec]=param[2];
1083       sec[isec]=isec;
1084       vecN[isec]=npoints;
1085       vecChi2[isec]=TMath::Sqrt(chi2/npoints);
1086
1087       TStatToolkit::FitPlane(tree,"0.264*CETmean.fElements", fstringG.Data(),Form("sector==%d",isec)+cutAll, chi2,npoints,param,covar,-1,0, 10000000, kFALSE);
1088       if (npoints<3) continue;
1089       printf("Sector=%d\n",isec);
1090       vecA0[isec]=param[0];
1091       vecALy[isec]=param[1];
1092       vecALx[isec]=param[2];
1093       vecAChi2[isec]=TMath::Sqrt(chi2/npoints);
1094     }
1095     isCalc=kTRUE;
1096     //
1097     TString  fstringRef="";              // global fit
1098     fstringRef+="gx.fElements++";
1099     fstringRef+="gy.fElements++";  
1100     fstringRef+="lx.fElements++";
1101     TStatToolkit::FitPlane(tree,"0.264*dt", fstringG.Data(),"(sector%36)<18"+cutAll, chi2,npoints,vecASide,covar,-1,0, 10000000, kFALSE);
1102     TStatToolkit::FitPlane(tree,"0.264*dt", fstringG.Data(),"(sector%36)>=18"+cutAll, chi2,npoints,vecCSide,covar,-1,0, 10000000, kFALSE);
1103
1104   }
1105   (*fPcstream)<<"dcs"<<     // CE information
1106     "CETSector.="  << &sec      << // sector numbers
1107     "CETRefA.="    << &vecASide << // diff to reference A side
1108     "CETRefC.="    << &vecCSide << // diff to reference C side
1109     //                      // fit in respect to reference data
1110     "CETRef0.="    << &vec0     << // offset change
1111     "CETRefY.="    << &vecLy    << // slope y change - rad
1112     "CETRefX.="    << &vecLx    << // slope x change - rad
1113     "CETRefChi2.=" << &vecChi2  << // chi2 (rms in cm)
1114     "CETRefN.="    << &vecN     << //number of accepted points
1115     //                       // fit in respect per mean per side
1116     "CET0.="       << &vecA0    << // offset change
1117     "CETY.="       << &vecALy   << // slope y change - rad
1118     "CETX.="       << &vecALx   << // slope x change - rad
1119     "CETChi2.="    << &vecAChi2  ; // chi2 (rms in cm)
1120 }
1121
1122 void AliTPCcalibSummary::ProcessPulserRef(){
1123   //
1124   // Get fit of residuals if Pulser in respect with reference
1125   // data
1126   //
1127   static TVectorD  sec(72);
1128   static TVectorD vec0(72);
1129   static TVectorD vecLy(72);
1130   static TVectorD vecLx(72);
1131   static TVectorD vecChi2(72);
1132   static TVectorD vecN(72);
1133   //
1134   static TVectorD vecA0(72);
1135   static TVectorD vecALy(72);
1136   static TVectorD vecALx(72);
1137   static TVectorD vecAChi2(72);
1138   static Bool_t isCalc=kFALSE;
1139   
1140   TFile f("calPads.root");
1141   TFile fref("calPadsRef.root");
1142   TTree * tree = (TTree*)f.Get("calPads");
1143   TTree * treeRef = (TTree*)fref.Get("calPads");
1144   tree->AddFriend(treeRef,"R");
1145   
1146   tree->SetAlias("inPulser","(abs(PulserTmean.fElements-PulserTmean_Median)<1.5&&abs(PulserTrms.fElements-PulserTrms_Median)<0.2)");  // outlyerTrms
1147   tree->SetAlias("inPulserR","(abs(R.PulserTmean.fElements-R.PulserTmean_Median)<1.5&&abs(R.PulserTrms.fElements-R.PulserTrms_Median)<0.2)");  // outlyerTrms
1148   //
1149   if (!isCalc){
1150     // make fits only once
1151     TStatToolkit toolkit;
1152     Double_t chi2=0;
1153     Int_t    npoints=0;
1154     TVectorD param;
1155     TMatrixD covar;
1156     tree->SetAlias("dt","PulserTmean.fElements-R.PulserTmean.fElements");
1157     TCut cutAll ="inPulser&&inPulserR"; 
1158     TString  fstringG="";              // global part
1159     fstringG+="ly.fElements++";
1160     fstringG+="(lx.fElements-134.)++";  
1161     for (Int_t isec=0; isec<72; isec++){
1162       TStatToolkit::FitPlane(tree,"dt", fstringG.Data(),Form("sector==%d",isec)+cutAll, chi2,npoints,param,covar,-1,0, 10000000, kFALSE);
1163       if (npoints<3) continue;
1164       printf("Setor=%d\n",isec);
1165       vec0[isec]=param[0];
1166       vecLy[isec]=param[1];
1167       vecLx[isec]=param[2];
1168       sec[isec]=isec;
1169       vecN[isec]=npoints;
1170
1171       TStatToolkit::FitPlane(tree,"PulserTmean.fElements", fstringG.Data(),Form("sector==%d",isec)+cutAll, chi2,npoints,param,covar,-1,0, 10000000, kFALSE);
1172       if (npoints<3) continue;
1173       printf("Setor=%d\n",isec);
1174       vecA0[isec]=param[0];
1175       vecALy[isec]=param[1];
1176       vecALx[isec]=param[2];
1177       vecAChi2[isec]=TMath::Sqrt(chi2/npoints);
1178     }
1179
1180     isCalc=kTRUE;
1181   }
1182   (*fPcstream)<<"dcs"<<     // Pulser information
1183     "PulserTSector.="  << &sec     << // sector numbers
1184     //                      // fit in respect to reference
1185     "PulserTRef0.="    << &vec0    << // offset change
1186     "PulserTRefY.="    << &vecLy   << // slope y change - rad
1187     "PulserTRefX.="    << &vecLx   << // slope x change - rad
1188     "PulserTRefChi2.=" << &vecChi2 << // chi2 (rms in cm)
1189     "PulserTRefN.="    << &vecN    << //number of accepted points
1190     //                       // fit in respect per mean per side
1191     "PulserT0.="       << &vecA0   << // offset change
1192     "PulserTY.="       << &vecALy  << // slope y change - rad
1193     "PulserTX.="       << &vecALx  << // slope x change - rad
1194     "PulserTChi2.="    << &vecAChi2 ; // chi2 (rms in cm)
1195 }
1196
1197 void AliTPCcalibSummary::ProcessCurrent(Int_t irun, Int_t itime){
1198   //
1199   // Dump current 
1200   //
1201   //variables to export
1202   //
1203   static TObjArray *currentArray=new TObjArray(72);   // current graphs
1204   static TObjArray *currentArray2=new TObjArray(72);  // current graphs to export
1205   //
1206   static TVectorD currentIROC(36);                    // current snapshots
1207   static TVectorD currentOROC(36); 
1208   static TVectorF sector(72);                         //
1209   static Double_t medcurIROC = 0;
1210   static Double_t medcurOROC = 0;
1211   //
1212   static TVectorF minROC(72);                         // current mean +-5 minutes
1213   static TVectorF maxROC(72);
1214   static TVectorF meanROC(72);
1215   static TVectorF medianROC(72);
1216   static Double_t meanIIROC=0;
1217   static Double_t meanIOROC=0;
1218   static Double_t medianIIROC=0;
1219   static Double_t medianIOROC=0;
1220   //
1221   AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(irun);
1222   //
1223   for(Int_t j=1; j<36; j++) currentIROC[j-1] = fCalibDB->GetChamberHighVoltage(irun, j,itime,-1,kTRUE);
1224   for(Int_t j=36; j<72; j++) currentOROC[j-36] = fCalibDB->GetChamberHighVoltage(irun, j,itime,-1,kTRUE);
1225   medcurIROC = TMath::Median(36, currentIROC.GetMatrixArray());
1226   medcurOROC = TMath::Median(36, currentOROC.GetMatrixArray());
1227
1228
1229   if (currentArray->At(0)==0){
1230     for (Int_t isec=0; isec<72; isec++){
1231       TString sensorName="";
1232       const char* sideName=(isec%36<18) ? "A":"C";
1233       if (isec<36){
1234         //IROC
1235         sensorName=Form("TPC_ANODE_I_%s%02d_IMEAS",sideName,isec%18);
1236       }else{
1237         //OROC
1238         sensorName=Form("TPC_ANODE_O_%s%02d_0_IMEAS",sideName,isec%18);
1239       }      
1240     
1241       AliDCSSensor *sensor = 0;
1242       if (voltageArray) sensor= voltageArray->GetSensor(sensorName);   
1243       TGraph *gr=0;
1244       if (!sensor) gr=new TGraph(1);
1245       else{
1246         if (!sensor->GetGraph()) gr=new TGraph(1);
1247         else{
1248           gr=sensor->GetGraph();
1249           Double_t startTime=sensor->GetStartTime();
1250           Double_t * time = new Double_t[gr->GetN()];
1251           for (Int_t ip=0; ip<gr->GetN(); ip++){ time[ip]= (gr->GetX()[ip]*3600.)+startTime;}     
1252           gr=new TGraph(gr->GetN(), time, gr->GetY());
1253           delete [] time;
1254         }      
1255       }
1256       gr->Sort();
1257       currentArray->AddAt(gr, isec);
1258       currentArray->AddAt(gr->Clone(), isec);
1259     }
1260   }
1261
1262
1263   for (Int_t isec=0; isec<72; isec++){
1264     sector[isec]=isec;
1265     TGraph * gr = (TGraph*)currentArray->At(isec);
1266     TGraph * graph2 = (TGraph*)currentArray2->At(isec);    
1267     Int_t firstBin= TMath::BinarySearch(gr->GetN(), gr->GetX(), itime-300.)-2;
1268     Int_t lastBin= TMath::BinarySearch(gr->GetN(), gr->GetX(), itime+300.)+2;
1269     if (firstBin<0) firstBin=0;
1270     if (lastBin>=gr->GetN()) lastBin=gr->GetN()-1;
1271     //
1272     if (firstBin<lastBin){
1273       //
1274       minROC[isec]=TMath::MinElement(lastBin-firstBin, &(gr->GetY()[firstBin]));
1275       maxROC[isec]=TMath::MaxElement(lastBin-firstBin, &(gr->GetY()[firstBin]));
1276       meanROC[isec]=TMath::Mean(lastBin-firstBin, &(gr->GetY()[firstBin]));
1277       medianROC[isec]=TMath::Median(lastBin-firstBin, &(gr->GetY()[firstBin]));
1278       graph2 = new TGraph(lastBin-firstBin, &(gr->GetX()[firstBin]), &(gr->GetY()[firstBin]));
1279       delete currentArray2->At(isec);
1280       currentArray2->AddAt(graph2,isec);
1281     }
1282     (*fPcstream)<<"dcs"<<     // current information
1283       Form("current%d.=",isec)<<graph2;
1284   }     
1285   meanIIROC=TMath::Mean(36, &(meanROC.GetMatrixArray()[0]));
1286   meanIOROC=TMath::Mean(36, &(meanROC.GetMatrixArray()[36]));
1287   medianIIROC=TMath::Median(36, &(meanROC.GetMatrixArray()[0]));
1288   medianIOROC=TMath::Median(36, &(meanROC.GetMatrixArray()[36]));
1289   //
1290   (*fPcstream)<<"dcs"<<     // current information
1291     "isec.="        << &sector      << // sector number
1292     "IIROC.="       << &currentIROC << // current sample at given moment
1293     "IOROC.="       << &currentOROC << // current sample at given moment
1294     "medianIIROC="  << medcurIROC   << // median at given moment
1295     "medianIOROC="  << medcurOROC   << // median at given moment
1296     //
1297     "minIROC.="     << &minROC      << // minimum in +-5 min
1298     "maxIROC.="     << &maxROC      << // maximum in +-5 min
1299     "meanIROC.="    << &meanROC     << // mean in +-5 min
1300     "medianIROC.="  << &medianROC   << // median in +-5 min
1301     "meanIIROC5="   << meanIIROC    << // mean current in IROC +-5 minutes
1302     "meanIOROC5="   << meanIOROC    << // mean current in OROC
1303     "medianIIROC5=" << medianIIROC  << // median current in IROC
1304     "medianIOROC5=" << medianIOROC   ; // medianan current in OROC
1305    
1306
1307   (*fPcstream)<<"current"<<     // current information
1308     "time="         << itime        <<
1309     "isec.="        << &sector      << // sector number
1310     "IIROC.="       << &currentIROC << // current sample at given moment
1311     "IOROC.="       << &currentOROC << // current sample at given moment
1312     "medianIIROC="  << medcurIROC   << // median at given moment
1313     "medianIOROC="  << medcurOROC   << // median at given moment
1314     //
1315     "minIROC.="     << &minROC      << // minimum in +-5 min
1316     "maxIROC.="     << &maxROC      << // maximum in +-5 min
1317     "meanIROC.="    << &meanROC     << // mean in +-5 min
1318     "medianIROC.="  << &medianROC   << // median in +-5 min
1319     "meanIIROC5="   << meanIIROC    << // mean current in IROC +-5 minutes
1320     "meanIOROC5="   << meanIOROC    << // mean current in OROC
1321     "medianIIROC5=" << medianIIROC  << // median current in IROC
1322     "medianIOROC5=" << medianIOROC  << // medianan current in OROC
1323     "\n";
1324
1325 }
1326
1327
1328
1329
1330 // TCanvas * DrawCEDiff(TTree * tree){
1331   
1332 //   TCanvas *canvasIO = new TCanvas("canvasCEIO","canvasCEIO");
1333 //   canvasIO->Divide(6,6);
1334 //   for (Int_t isec=0; isec<36; isec++){
1335 //     canvasIO->cd(isec+1);
1336 //     dcs->Draw(Form("CET0.fElements[%d]-CET0.fElements[%d]",isec+36,isec),Form("abs(CETRef0.fElements[%d])<0.3",isec),"");
1337 //     printf("%d\t%f\t%f\n",isec,dcs->GetHistogram()->GetMean(),dcs->GetHistogram()->GetRMS());
1338 //   }
1339
1340 // }