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