]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCcalibSummary.cxx
Forgotten commit
[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 //
a3b590cf 176 Int_t dtime = TMath::Max((endTime-startTime)/20,10);
12f69174 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 //
abb20887 232 //QA data processing
233 //
234 TVectorD vQaOcc;
235 TVectorD vQaQtot;
236 TVectorD vQaQmax;
237 fDButil->ProcessQAData(vQaOcc, vQaQtot, vQaQmax);
238 //
12f69174 239 //calibration Pulser data processing
240 //
241 Int_t nOffChannels=0;
242 TVectorD vTimePulser;
243 nOffChannels=fDButil->GetNPulserOutliers();
244 fDButil->ProcessPulser(vTimePulser);
245 //
246 //ALTRO data
247 //
248 Int_t nMasked=0;
249 fDButil->ProcessALTROConfig(nMasked);
250 //
251 //Calib RAW data
252 //
253 Int_t nFailL1=-1;
254 if (fCalibDB->GetCalibRaw()) nFailL1=fCalibDB->GetCalibRaw()->GetNFailL1Phase();
255 //
256 //production information
257 //
258 Int_t nalien=0,nRawAlien=0,nlocal=0,nRawLocal=0;
259 //run type
260 TObjString runType(AliTPCcalibDB::GetRunType(irun).Data());
261 //
262 //
263 //
264
265 for (Int_t itime=startTime; itime<endTime; itime+=dtime){
266 //
267 TTimeStamp tstamp(itime);
268 Float_t valuePressure = fCalibDB->GetPressure(tstamp,irun,0);
269 Float_t valuePressure2 = fCalibDB->GetPressure(tstamp,irun,1);
270 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative((UInt_t)itime,irun,0);
271 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative((UInt_t)itime,irun,1);
272 //temperature fits
273 TLinearFitter * fitter = 0;
274 TVectorD vecTemp[10];
275 for (Int_t itype=0; itype<5; itype++)
276 for (Int_t iside=0; iside<2; iside++){
277 fitter= tempMap->GetLinearFitter(itype,iside,tstamp);
278 if (!fitter) continue;
279 fitter->Eval();
280 fitter->GetParameters(vecTemp[itype+iside*5]);
281 delete fitter;
282 }
283 //
284 //measured skirt temperatures
285 //
286 TVectorD vecSkirtTempA(18);
287 TVectorD vecSkirtTempC(18);
288 Int_t nsenTemp=tempArray->NumSensors();
289 for (Int_t isenTemp=0;isenTemp<nsenTemp;++isenTemp){
290 AliTPCSensorTemp *senTemp=(AliTPCSensorTemp*)tempArray->GetSensorNum(isenTemp);
291 if (senTemp->GetType()!=3) continue;
292 if (TMath::Sqrt(senTemp->GetX()*senTemp->GetX()+senTemp->GetY()*senTemp->GetY())<100) continue; //only skirt, outer FC vessel
293 Double_t val=senTemp->GetValue(tstamp);
294 if (senTemp->GetSide()==0)
295 vecSkirtTempA[senTemp->GetSector()]=val;
296 else
297 vecSkirtTempC[senTemp->GetSector()]=val;
298 }
299 //
300 //goofie data
301 //
302 TVectorD vecGoofie;
303 if (goofieArray){
304 vecGoofie.ResizeTo(goofieArray->NumSensors());
305 for (Int_t isensor=0; isensor<goofieArray->NumSensors();isensor++){
306 AliDCSSensor *gsensor = goofieArray->GetSensor(isensor);
307 if (gsensor){
308 vecGoofie[isensor] = gsensor->GetValue(tstamp);
309 }
310 }
311 } else {
312 vecGoofie.ResizeTo(19);
313 }
314 //
315 TVectorD voltagesIROC(36);
a3b590cf 316 TVectorD voltagesOROC(36);
12f69174 317 for(Int_t j=1; j<36; j++) voltagesIROC[j-1] = fCalibDB->GetChamberHighVoltage(irun, j,itime);
318 for(Int_t j=36; j<72; j++) voltagesOROC[j-36] = fCalibDB->GetChamberHighVoltage(irun, j,itime);
319 Double_t voltIROC = TMath::Median(36, voltagesIROC.GetMatrixArray());
320 Double_t voltOROC = TMath::Median(36, voltagesOROC.GetMatrixArray());
321 //
322 Float_t coverIA=AliTPCcalibDB::GetCoverVoltage(irun,0,itime);
323 Float_t coverIC=AliTPCcalibDB::GetCoverVoltage(irun,18,itime);
324 Float_t coverOA=AliTPCcalibDB::GetCoverVoltage(irun,36,itime);
325 Float_t coverOC=AliTPCcalibDB::GetCoverVoltage(irun,54,itime);
326 Float_t skirtA=AliTPCcalibDB::GetSkirtVoltage(irun,0,itime);
327 Float_t skirtC=AliTPCcalibDB::GetSkirtVoltage(irun,18,itime);
328 Float_t ggOffA=AliTPCcalibDB::GetGGoffsetVoltage(irun,0,itime);
329 Float_t ggOffC=AliTPCcalibDB::GetGGoffsetVoltage(irun,18,itime);
330 //drift velocity
331 Float_t dvCorr=-5;
332 if (fitVdrift) dvCorr=fitVdrift->Eval(itime);
333
334 //tempMap->GetLinearFitter(0,0,itime);
335 (*fPcstream)<<"dcs"<<
336 "run="<<irun<<
337 "time="<<itime<<
338 "startTimeGRP="<<startTimeGRP<<
339 "stopTimeGRP="<<stopTimeGRP<<
340 //run type
341 "runType.="<<&runType<<
342 // voltage setting
343 "VIROC.="<<&voltagesIROC<<
344 "VOROC.="<<&voltagesOROC<<
345 "medianVIROC="<<voltIROC<<
346 "medianVOROC="<<voltOROC<<
347 "coverIA=" << coverIA <<
348 "coverIC=" << coverIC <<
349 "coverOA=" << coverOA <<
350 "coverOC=" << coverOC <<
351 "skirtA=" << skirtA <<
352 "skirtC=" << skirtC <<
353 "ggOffA=" << ggOffA <<
354 "ggOffC=" << ggOffC <<
355 //
356 "ptrel0="<<ptrelative0<< // deltaTP/TP - A side
357 "ptrel1="<<ptrelative1<< // deltaTP/TPC - C side
358 "goofie.="<<&vecGoofie<<
359 "goofieE.="<<&vecEntries<<
360 "goofieMean.="<<&vecMean<<
361 "goofieMedian.="<<&vecMedian<<
362 "goofieRMS.="<<&vecRMS<<
363 //
364 "press="<<valuePressure<<
365 "press2="<<valuePressure2<<
366 "temp00.="<<&vecTemp[0]<<
367 "temp10.="<<&vecTemp[1]<<
368 "temp20.="<<&vecTemp[2]<<
369 "temp30.="<<&vecTemp[3]<<
370 "temp40.="<<&vecTemp[4]<<
371 "temp01.="<<&vecTemp[5]<<
372 "temp11.="<<&vecTemp[6]<<
373 "temp21.="<<&vecTemp[7]<<
374 "temp31.="<<&vecTemp[8]<<
375 "temp41.="<<&vecTemp[9]<<
376 "tempSkirtA.="<<&vecSkirtTempA<<
377 "tempSkirtC.="<<&vecSkirtTempC;
378
379 ProcessDrift(irun, itime);
380 ProcessDriftCE(irun,itime);
381 ProcessDriftAll(irun,itime);
5f547d58 382 // ProcessKryptonTime(irun,itime);
12f69174 383 ProcessCTP(irun,itime);
384 ProcessAlign(irun,itime);
385 ProcessGain(irun,itime);
a3b590cf 386 //ProcessDriftCERef();
387 //ProcessPulserRef();
388 ProcessCurrent(irun,itime);
12f69174 389
390
391 (*fPcstream)<<"dcs"<<
392 //noise data
393 "meanNoise.="<<&vNoiseMean<<
394 "meanNoiseSen.="<<&vNoiseMeanSenRegions<<
395 "rmsNoise.="<<&vNoiseRMS<<
396 "rmsNoiseSen.="<<&vNoiseRMSSenRegions<<
397 "zeroNoise="<<nonMaskedZero<<
abb20887 398 "nNaN="<<nNaN<<
399 //QA data
400 "occQA.=" << &vQaOcc <<
401 "qQA.=" << &vQaQtot <<
402 "qmaxQA.=" << &vQaQmax <<
12f69174 403 //pulser data
404 "timePulser.=" << &vTimePulser <<
405 "nOffPulser="<<nOffChannels<<
406 //altro data
407 "nMasked="<<nMasked<<
408 //ce data -Jens version
409 "CEfitA.="<<&fitResultsA<<
410 "CEfitC.="<<&fitResultsC<<
411 "nmaskedCE="<<nmaskedCE<<
412 "chi2ACE="<<chi2ACE<<
413 "chi2CCE="<<chi2CCE<<
414 //ce data new - MI version
415 "CEfitAMI.="<<&fitCEResultsA<<
416 "CEfitCMI.="<<&fitCEResultsC<<
417 "chi2CEA="<<chi2CEA<<
418 "chi2CEC="<<chi2CEC<<
419 //
420 //ce graph data
421 "CEgrTEntries.="<<&vecTEntries<<
422 "CEgrTMean.="<<&vecTMean<<
423 "CEgrTRMS.="<<&vecTRMS<<
424 "CEgrTMedian.="<<&vecTMedian<<
425 "CEgrQEntries.="<<&vecQEntries<<
426 "CEgrQMean.="<<&vecQMean<<
427 "CEgrQRMS.="<<&vecQRMS<<
428 "CEgrQMedian.="<<&vecQMedian<<
429 "CEgrDriftA="<<driftTimeA<<
430 "CEgrDriftC="<<driftTimeC<<
431 //calib raw data
432 "nFailL1="<<nFailL1<<
433 // b field
434 "Bz="<< bz <<
435 "L3polarity="<<l3pol<<
436 // production information
437 "nalien="<<nalien<<
438 "nRawAlien="<<nRawAlien<<
439 "nlocal="<<nlocal<<
440 "nRawLocal="<<nRawLocal<<
441 //comparisons with ref data
442 "pedestalDeviations.="<<&pedestalDeviations<<
443 "noiseDeviations.="<<&noiseDeviations<<
444 "pulserQdeviations.="<<&pulserQdeviations<<
445 // "pulserVarQMean="<<varQMean<<
446 "pulserNpadsOutOneTB="<<npadsOutOneTB<<
447 "pulserNpadsOffAdd="<<npadsOffAdd<<
448 "driftCorrCosmAll="<<dvCorr<<
449 "\n";
450 }//end run loop
451}
452
453
454
455
456
457
458void AliTPCcalibSummary::ProcessDrift(Int_t run, Int_t timeStamp){
459 //
460 // dump drift calibration data to the tree
461 //
462 TObjArray *array =fCalibDB->GetTimeVdriftSplineRun(run);
463 TGraphErrors *laserA[3]={0,0,0};
464 TGraphErrors *laserC[3]={0,0,0};
465 TGraphErrors *cosmicAll=0;
466 static Double_t vlaserA[3]={0,0,0};
467 static Double_t vlaserC[3]={0,0,0};
468 static Double_t vcosmicAll=0;
469 static Double_t vdrift1=0;
470 vdrift1=fCalibDB->GetVDriftCorrectionTime(timeStamp,run,0,1);
471
472 if (array){
473 laserA[0]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DELAY_LASER_ALL_A");
474 laserA[1]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_A");
475 laserA[2]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_A");
476 laserC[0]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DELAY_LASER_ALL_C");
477 laserC[1]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_DRIFT_LASER_ALL_C");
478 laserC[2]=(TGraphErrors*)array->FindObject("GRAPH_MEAN_GLOBALYGRADIENT_LASER_ALL_C");
479 cosmicAll =(TGraphErrors*)array->FindObject("TGRAPHERRORS_MEAN_VDRIFT_COSMICS_ALL");
480 }
481 if (laserA[0]) vlaserA[0]= AliTPCcalibDButil::EvalGraphConst(laserA[0],timeStamp);
482 if (laserA[1]) vlaserA[1]= AliTPCcalibDButil::EvalGraphConst(laserA[1],timeStamp);
483 if (laserA[2]) vlaserA[2]= AliTPCcalibDButil::EvalGraphConst(laserA[2],timeStamp);
484 if (laserC[0]) vlaserC[0]= AliTPCcalibDButil::EvalGraphConst(laserC[0],timeStamp);
485 if (laserC[1]) vlaserC[1]= AliTPCcalibDButil::EvalGraphConst(laserC[1],timeStamp);
486 if (laserC[2]) vlaserC[2]= AliTPCcalibDButil::EvalGraphConst(laserC[2],timeStamp);
487 if (cosmicAll) vcosmicAll= AliTPCcalibDButil::EvalGraphConst(cosmicAll,timeStamp);
488 (*fPcstream)<<"dcs"<<
489 "vlaserA0="<<vlaserA[0]<<
490 "vlaserA1="<<vlaserA[1]<<
491 "vlaserA2="<<vlaserA[2]<<
492 "vlaserC0="<<vlaserC[0]<<
493 "vlaserC1="<<vlaserC[1]<<
494 "vlaserC2="<<vlaserC[2]<<
495 "vcosmicAll="<<vcosmicAll<<
496 "vdrift1="<<vdrift1;
497
498 //
499 // define distance to measurement
500 //
501 static Double_t dlaserA=0;
502 static Double_t dlaserC=0;
503 static Double_t dcosmic=0;
504 static Double_t slaserA=0;
505 static Double_t slaserC=0;
506 static Double_t scosmic=0;
507 static Double_t vclaserA[3]={0,0,0};
508 static Double_t vclaserC[3]={0,0,0};
509 static Double_t vccosmicAll=0;
510 for (Int_t i=0;i<3;i++){
511 if (laserA[i]) AliTPCcalibDButil::GetNearest(laserA[i],timeStamp,dlaserA,vclaserA[i]);
512 if (laserC[i]) AliTPCcalibDButil::GetNearest(laserC[i],timeStamp,dlaserC,vclaserC[i]);
513 }
514 if (cosmicAll) AliTPCcalibDButil::GetNearest(cosmicAll,timeStamp,dcosmic,vccosmicAll);
515 (*fPcstream)<<"dcs"<<
516 "vclaserA0="<<vclaserA[0]<<
517 "vclaserA1="<<vclaserA[1]<<
518 "vclaserA2="<<vclaserA[2]<<
519 "vclaserC0="<<vclaserC[0]<<
520 "vclaserC1="<<vclaserC[1]<<
521 "vclaserC2="<<vclaserC[2]<<
522 "vccosmicAll="<<vccosmicAll<<
523 "dlaserA="<<dlaserA<<
524 "dlaserC="<<dlaserC<<
525 "dcosmic="<<dcosmic<<
526 "slaserA="<<slaserA<<
527 "slaserC="<<slaserC<<
528 "scosmic="<<scosmic;
529}
530
531void AliTPCcalibSummary::ProcessDriftCE(Int_t run,Int_t timeStamp){
532 //
533 // dump drift calibration data CE
534 //
535 TObjArray *arrT=fCalibDB->GetCErocTtime();
536 AliTPCParam *param=fCalibDB->GetParameters();
537 static TVectorD tdriftCE(74);
538 static TVectorD tndriftCE(74);
539 static TVectorD vdriftCE(74);
540 static TVectorD tcdriftCE(74);
541 static TVectorD tddriftCE(74);
542 static Double_t ltime0A;
543 static Double_t ltime0C;
544 //
545 //
546 //
547 ltime0A = fDButil->GetLaserTime0(run,timeStamp,36000,0);
548 ltime0C = fDButil->GetLaserTime0(run,timeStamp,36000,1);
549 //
550 for (Int_t i=0; i<arrT->GetEntries();i++){
551 tdriftCE[i]=0;
552 vdriftCE[i]=0;
553 TGraph *graph = (TGraph*)arrT->At(i);
554 if (!graph) continue;
555 tdriftCE[i]=AliTPCcalibDButil::EvalGraphConst(graph,timeStamp);
556 Double_t deltaT,gry;
557 AliTPCcalibDButil::GetNearest(graph,timeStamp,deltaT,gry);
558 tndriftCE[i]=graph->GetN();
559 tcdriftCE[i]=gry;
560 tddriftCE[i]=deltaT;
561 if (i%36<18){
562 vdriftCE[i] =(param->GetZLength(i)/(tdriftCE[i]*param->GetTSample()*(1.-ltime0A)-param->GetL1Delay()))/param->GetDriftV();
563 }else{
564 vdriftCE[i] =(param->GetZLength(i)/(tdriftCE[i]*param->GetTSample()*(1.-ltime0A)-param->GetL1Delay()))/param->GetDriftV();
565 }
566 }
567 // export values
568 (*fPcstream)<<"dcs"<<
569 "tdriftCE.="<<&tdriftCE<< // arrival time
570 "vdriftCE.="<<&vdriftCE<< // derived drift velocity per chamber
571 "tndriftCE.="<<&tndriftCE<< // number of points (chambers)
572 "tcdriftCE.="<<&tcdriftCE<< // constant evaluation - nearest point used
573 "tddriftCE.="<<&tddriftCE<< // distance to closest measuement
574 "ltime0A="<<ltime0A<< // laser offset expected in reconstruction
575 "ltime0C="<<ltime0C; // laser offset expected in reconstruction
576 }
577
578
579void AliTPCcalibSummary::ProcessDriftAll(Int_t run,Int_t timeStamp){
580 //
581 // dump drift calibration data all calibrations form DB util
582 // test of utils
583 static Double_t vdriftCEA=0, vdriftCEC=0, vdriftCEM=0;
584 static Double_t vdriftLTA=0, vdriftLTC=0, vdriftLTM=0;
78f17711 585 static Double_t vdriftLTAon=0, vdriftLTCon=0, vdriftLTMon=0;
12f69174 586 static Double_t vdriftITS=0;
587 static Double_t vdriftP=0;
78f17711 588 static Double_t dcea=0, dcec=0, dcem=0, dla=0,dlc=0,dlm=0, dlaon=0,dlcon=0,dlmon=0, dp=0;
12f69174 589 static Double_t dits=0;
590 static Double_t ltime0A;
591 static Double_t ltime0C;
592 static Double_t ctime0;
593 static Double_t vdrift1=0;
594 vdrift1= fCalibDB->GetVDriftCorrectionTime(timeStamp,run,0,1);
595 vdriftP = fDButil->GetVDriftTPC(dp, run, timeStamp, 86400, 3600,0);
596 ctime0 = AliTPCcalibDButil::GetTriggerOffsetTPC(run,timeStamp, 36000, 3600,0);
597 //
598 vdriftCEA= fDButil->GetVDriftTPCCE(dcea,run,timeStamp,36000,0);
599 vdriftCEC= fDButil->GetVDriftTPCCE(dcec,run,timeStamp,36000,1);
600 vdriftCEM= fDButil->GetVDriftTPCCE(dcem,run,timeStamp,36000,2);
601 //
602 vdriftLTA= fDButil->GetVDriftTPCLaserTracks(dla,run,timeStamp,36000,0);
603 vdriftLTC= fDButil->GetVDriftTPCLaserTracks(dlc,run,timeStamp,36000,1);
604 vdriftLTM= fDButil->GetVDriftTPCLaserTracks(dlm,run,timeStamp,36000,2);
605 //
78f17711 606 vdriftLTAon= fDButil->GetVDriftTPCLaserTracksOnline(dlaon,run,timeStamp,36000,0);
607 vdriftLTCon= fDButil->GetVDriftTPCLaserTracksOnline(dlcon,run,timeStamp,36000,1);
608 vdriftLTMon= fDButil->GetVDriftTPCLaserTracksOnline(dlmon,run,timeStamp,36000,2);
609 //
12f69174 610 vdriftITS= fDButil->GetVDriftTPCITS(dits, run,timeStamp);
611 //
612 ltime0A = fDButil->GetLaserTime0(run,timeStamp,36000,0);
613 ltime0C = fDButil->GetLaserTime0(run,timeStamp,36000,1);
614
615 (*fPcstream)<<"dcs"<<
616 //
617 "vdriftCEA="<<vdriftCEA<< // drift velocity CE
618 "vdriftCEC="<<vdriftCEC<<
619 "vdriftCEM="<<vdriftCEM<<
620 "dcea="<<dcea<<
621 "dcec="<<dcec<<
622 "dcem="<<dcem<<
623 "vdriftLTA="<<vdriftLTA<< // drift velocity laser tracks
624 "vdriftLTC="<<vdriftLTC<<
625 "vdriftLTM="<<vdriftLTM<<
626 "dla="<<dla<<
627 "dlc="<<dlc<<
628 "dlm="<<dlm<<
78f17711 629 "vdriftLTAon="<<vdriftLTAon<< // drift velocity laser tracks and CE from online algorithm
630 "vdriftLTCon="<<vdriftLTCon<<
631 "vdriftLTMon="<<vdriftLTMon<<
632 "dlaOn="<<dlaon<<
633 "dlcOn="<<dlcon<<
634 "dlmOn="<<dlmon<<
12f69174 635 //
636 //
637 "vdriftITS="<<vdriftITS<<
638 "dits="<<dits<<
639 "ctime0="<<ctime0<<
640 "vdriftP="<<vdriftP<< // drift velocity comsic
641 "dp="<<dp<<
642 "vdrift1="<<vdrift1; // combined drift velocity
643
644}
645
646
647
648void AliTPCcalibSummary::ProcessKryptonTime(Int_t run, Int_t timeStamp){
649 //
650 // Dumping krypton calibration results
651 //
652 static TObjArray * krArray=0;
653 if (!krArray) {
654 AliCDBEntry* entry = AliCDBManager::Instance()->Get("TPC/Calib/TimeGainKrypton", run);
655 if (entry){
656 krArray = (TObjArray*)entry->GetObject();
657 }
658 }
659 static TVectorD krMean(74);
660 static TVectorD krErr(74);
661 static TVectorD krDist(74);
662 TGraphErrors *gr=0;
663 Double_t deltaT=0,gry=0;
664 if (krArray){
665 for (Int_t isec=0; isec<72; isec++){
666 krMean[isec]=0;
667 krDist[isec]=0;
668 krErr[isec]=0;
669 gr=(TGraphErrors*)krArray->At(isec);
670 if (gr) {
5647625c 671 krMean[isec]=AliTPCcalibDButil::EvalGraphConst(gr,timeStamp);
12f69174 672 AliTPCcalibDButil::GetNearest(gr, timeStamp,deltaT,gry);
673 krDist[isec]=deltaT;
674 }
675 if (72+isec<krArray->GetEntries()) {
676 gr=(TGraphErrors*)krArray->At(72+isec);
5647625c 677 if (gr) krErr[isec]=AliTPCcalibDButil::EvalGraphConst(gr,timeStamp);
12f69174 678 }
679 }
680 krMean[72]= TMath::Median(36,krMean.GetMatrixArray());
681 krMean[73]= TMath::Median(36,&(krMean.GetMatrixArray()[36]));
682 krErr[72]= TMath::Median(36,krErr.GetMatrixArray());
683 krErr[73]= TMath::Median(36,&(krErr.GetMatrixArray()[36]));
684 }
685 (*fPcstream)<<"dcs"<<
686 "krMean.="<<&krMean<<
687 "krErr.="<<&krErr<<
688 "krDist.="<<&krDist;
689}
690
691void AliTPCcalibSummary::ProcessCTP(Int_t irun, Int_t timeStamp){
692 //
693 //
694 //
695 static TClonesArray *pcarray = new TClonesArray("AliCTPInputTimeParams",1);
696 static AliCTPTimeParams *ctpParams =0;
697 ctpParams = fCalibDB->GetCTPTimeParams(); //
698 const TObjArray *arr = ctpParams->GetInputTimeParams();
699 pcarray->ExpandCreateFast(TMath::Max(arr->GetEntries(),1));
700 for (Int_t i=0; i<arr->GetEntries(); i++){
701 new ((*pcarray)[i]) AliCTPInputTimeParams(*((AliCTPInputTimeParams*)arr->At(i)));
702 }
703 (*fPcstream)<<"ctp"<<
704 "run="<<irun<<
705 "time="<<timeStamp<<
706 "ctpP.="<<ctpParams<<
707 "ctpArr="<<pcarray<<
708 "\n";
709 (*fPcstream)<<"dcs"<<
710 "ctpP.="<<ctpParams<<
711 "ctpArr="<<pcarray;
712}
713
714void AliTPCcalibSummary::ProcessGain(Int_t irun, Int_t timeStamp){
715 //
716 // Dump gain related information to the tree
717 //
718 static Float_t gainCosmic = 0;
719 static Float_t gainMIP = 0;
5647625c 720 static Float_t attachMIP = 0;
721 static Double_t dMIP=0;
722 Double_t dummy=0;
12f69174 723 TObjArray * gainSplines = fCalibDB->GetTimeGainSplinesRun(irun);
724 if (gainSplines) {
725 TGraphErrors * graphMIP = (TGraphErrors *) gainSplines->FindObject("TGRAPHERRORS_MEAN_GAIN_BEAM_ALL");
726 TGraphErrors * graphCosmic = (TGraphErrors *) gainSplines->FindObject("TGRAPHERRORS_MEAN_GAIN_COSMIC_ALL");
5647625c 727 TGraphErrors * graphAttach = (TGraphErrors *) gainSplines->FindObject("TGRAPHERRORS_MEAN_ATTACHMENT_BEAM_ALL");
12f69174 728 if (graphMIP) gainMIP = AliTPCcalibDButil::EvalGraphConst(graphMIP,timeStamp);
729 if (graphCosmic) gainCosmic = AliTPCcalibDButil::EvalGraphConst(graphCosmic,timeStamp);
5647625c 730 if (graphAttach) attachMIP = AliTPCcalibDButil::EvalGraphConst(graphAttach,timeStamp);
731 if (graphMIP) AliTPCcalibDButil::GetNearest(graphMIP, timeStamp, dMIP,dummy);
12f69174 732 }
733 // time dependence of gain
734 (*fPcstream)<<"dcs"<<
735 "gainMIP="<<gainMIP<<
5647625c 736 "attachMIP="<<attachMIP<<
737 "dMIP="<<dMIP<<
12f69174 738 "gainCosmic="<<gainCosmic;
739}
740
741
742void AliTPCcalibSummary::ProcessAlign(Int_t run, Int_t timeStamp){
743 //
744 // Proccess alignment
745 //
746 TString grnames[12]={"ALIGN_ITS", "ALIGN_ITSP", "ALIGN_ITSM", "ALIGN_ITSB",
747 "ALIGN_TRD", "ALIGN_TRDP", "ALIGN_TRDM","ALIGN_TRDB",
748 "ALIGN_TOF", "ALIGN_TOFP", "ALIGN_TOFM","ALIGN_TOFB"};
749 TString grpar[9]={"DELTAPSI", "DELTATHETA", "DELTAPHI",
750 "DELTAX", "DELTAY", "DELTAZ",
751 "DRIFTVD", "T0", "VDGY"};
752 static Double_t values[12*9];
753 static Double_t errs[12*9];
754
755 TObjArray *array =fCalibDB->GetTimeVdriftSplineRun(run);
756 TGraphErrors *gr=0;
757 for (Int_t idet=0; idet<12; idet++){
758 for (Int_t ipar=0; ipar<9; ipar++){
759 TString grName=grnames[idet];
760 grName+="_TPC_";
761 grName+=grpar[ipar];
762 if (array){
763 gr = (TGraphErrors*)array->FindObject(grName.Data());
764 }
765 values[9*idet+ipar]=0;
766 errs[9*idet+ipar]=0;
5647625c 767 if (gr) values[9*idet+ipar]=AliTPCcalibDButil::EvalGraphConst(gr,timeStamp);
12f69174 768 (*fPcstream)<<"dcs"<<
769 Form("%s=",grName.Data())<<values[9*idet+ipar];
770 (*fPcstream)<<"align"<<
771 Form("%s=",grName.Data())<<values[9*idet+ipar];
772 }
773 }
774 (*fPcstream)<<"align"<<
775 "time="<<timeStamp<<
776 "run="<<run<<
777 "\n";
778}
779
780
af6a50bb 781void AliTPCcalibSummary::ProcessDriftCERef(){
782 //
783 // Get fit of residuals if CE in respect with reference
784 // data
785 //
786 static TVectorD sec(72);
787 static TVectorD vec0(72);
788 static TVectorD vecLy(72);
789 static TVectorD vecLx(72);
790 static TVectorD vecChi2(72);
791 static TVectorD vecN(72);
792 //
793 static TVectorD vecA0(72);
794 static TVectorD vecALy(72);
795 static TVectorD vecALx(72);
796 static TVectorD vecAChi2(72);
797 //
5647625c 798 static TVectorD vecASide(4);
799 static TVectorD vecCSide(4);
af6a50bb 800 static Bool_t isCalc=kFALSE;
801
802 TFile f("calPads.root");
803 TFile fref("calPadsRef.root");
804 TTree * tree = (TTree*)f.Get("calPads");
805 TTree * treeRef = (TTree*)fref.Get("calPads");
806 tree->AddFriend(treeRef,"R");
807 tree->SetAlias("inCE","((CEQmean.fElements>35)&&abs(CETmean.fElements)<1.5&&abs(CETrms.fElements/1.2-1)<0.2)"); // outlyerTrms
808 tree->SetAlias("inCER","((R.CEQmean.fElements>35)&&abs(R.CETmean.fElements)<1.5&&abs(R.CETrms.fElements/1.2-1)<0.2)"); // outlyerTrms
809 //
810 if (!isCalc){
811 // make fits only once
812 TStatToolkit toolkit;
813 Double_t chi2=0;
814 Int_t npoints=0;
815 TVectorD param;
816 TMatrixD covar;
817 tree->SetAlias("dt","CETmean.fElements-R.CETmean.fElements");
818 TCut cutAll ="inCE&&inCER&&abs(CETmean.fElements-R.CETmean.fElements)<0.5";
819 TString fstringG=""; // global part
820 fstringG+="ly.fElements++";
821 fstringG+="(lx.fElements-134.)++";
822 for (Int_t isec=0; isec<72; isec++){
5647625c 823 TStatToolkit::FitPlane(tree,"0.264*dt", fstringG.Data(),Form("sector==%d",isec)+cutAll, chi2,npoints,param,covar,-1,0, 10000000, kFALSE);
af6a50bb 824 if (npoints<3) continue;
825 printf("Sector=%d\n",isec);
826 vec0[isec]=param[0];
827 vecLy[isec]=param[1];
828 vecLx[isec]=param[2];
829 sec[isec]=isec;
830 vecN[isec]=npoints;
5647625c 831 vecChi2[isec]=TMath::Sqrt(chi2/npoints);
af6a50bb 832
5647625c 833 TStatToolkit::FitPlane(tree,"0.264*CETmean.fElements", fstringG.Data(),Form("sector==%d",isec)+cutAll, chi2,npoints,param,covar,-1,0, 10000000, kFALSE);
af6a50bb 834 if (npoints<3) continue;
835 printf("Sector=%d\n",isec);
836 vecA0[isec]=param[0];
837 vecALy[isec]=param[1];
838 vecALx[isec]=param[2];
839 vecAChi2[isec]=TMath::Sqrt(chi2/npoints);
af6a50bb 840 }
841 isCalc=kTRUE;
5647625c 842 //
843 TString fstringRef=""; // global fit
844 fstringRef+="gx.fElements++";
845 fstringRef+="gy.fElements++";
846 fstringRef+="lx.fElements++";
847 TStatToolkit::FitPlane(tree,"0.264*dt", fstringG.Data(),"(sector%36)<18"+cutAll, chi2,npoints,vecASide,covar,-1,0, 10000000, kFALSE);
848 TStatToolkit::FitPlane(tree,"0.264*dt", fstringG.Data(),"(sector%36)>=18"+cutAll, chi2,npoints,vecCSide,covar,-1,0, 10000000, kFALSE);
849
af6a50bb 850 }
851 (*fPcstream)<<"dcs"<< // CE information
852 "CETSector.="<<&sec<< // sector numbers
5647625c 853 "CETRefA.="<<&vecASide<< // diff to reference A side
854 "CETRefC.="<<&vecCSide<< // diff to reference C side
855 // // fit in respect to reference data
af6a50bb 856 "CETRef0.="<<&vec0<< // offset change
857 "CETRefY.="<<&vecLy<< // slope y change - rad
858 "CETRefX.="<<&vecLx<< // slope x change - rad
859 "CETRefChi2.="<<&vecChi2<< // chi2 (rms in cm)
860 "CETRefN.="<<&vecN<< //number of accepted points
861 // // fit in respect per mean per side
862 "CET0.="<<&vecA0<< // offset change
863 "CETY.="<<&vecALy<< // slope y change - rad
864 "CETX.="<<&vecALx<< // slope x change - rad
865 "CETChi2.="<<&vecAChi2; // chi2 (rms in cm)
866}
12f69174 867
af6a50bb 868void AliTPCcalibSummary::ProcessPulserRef(){
869 //
870 // Get fit of residuals if Pulser in respect with reference
871 // data
872 //
873 static TVectorD sec(72);
874 static TVectorD vec0(72);
875 static TVectorD vecLy(72);
876 static TVectorD vecLx(72);
877 static TVectorD vecChi2(72);
878 static TVectorD vecN(72);
879 //
880 static TVectorD vecA0(72);
881 static TVectorD vecALy(72);
882 static TVectorD vecALx(72);
883 static TVectorD vecAChi2(72);
884 static Bool_t isCalc=kFALSE;
885
886 TFile f("calPads.root");
887 TFile fref("calPadsRef.root");
888 TTree * tree = (TTree*)f.Get("calPads");
889 TTree * treeRef = (TTree*)fref.Get("calPads");
890 tree->AddFriend(treeRef,"R");
891
892 tree->SetAlias("inPulser","(abs(PulserTmean.fElements-PulserTmean_Median)<1.5&&abs(PulserTrms.fElements-PulserTrms_Median)<0.2)"); // outlyerTrms
893 tree->SetAlias("inPulserR","(abs(R.PulserTmean.fElements-R.PulserTmean_Median)<1.5&&abs(R.PulserTrms.fElements-R.PulserTrms_Median)<0.2)"); // outlyerTrms
894 //
895 if (!isCalc){
896 // make fits only once
897 TStatToolkit toolkit;
898 Double_t chi2=0;
899 Int_t npoints=0;
900 TVectorD param;
901 TMatrixD covar;
902 tree->SetAlias("dt","PulserTmean.fElements-R.PulserTmean.fElements");
903 TCut cutAll ="inPulser&&inPulserR";
904 TString fstringG=""; // global part
905 fstringG+="ly.fElements++";
906 fstringG+="(lx.fElements-134.)++";
907 for (Int_t isec=0; isec<72; isec++){
908 TStatToolkit::FitPlane(tree,"dt", fstringG.Data(),Form("sector==%d",isec)+cutAll, chi2,npoints,param,covar,-1,0, 10000000, kFALSE);
909 if (npoints<3) continue;
910 printf("Setor=%d\n",isec);
911 vec0[isec]=param[0];
912 vecLy[isec]=param[1];
913 vecLx[isec]=param[2];
914 sec[isec]=isec;
915 vecN[isec]=npoints;
916
917 TStatToolkit::FitPlane(tree,"PulserTmean.fElements", fstringG.Data(),Form("sector==%d",isec)+cutAll, chi2,npoints,param,covar,-1,0, 10000000, kFALSE);
918 if (npoints<3) continue;
919 printf("Setor=%d\n",isec);
920 vecA0[isec]=param[0];
921 vecALy[isec]=param[1];
922 vecALx[isec]=param[2];
923 vecAChi2[isec]=TMath::Sqrt(chi2/npoints);
924 }
5647625c 925
af6a50bb 926 isCalc=kTRUE;
927 }
928 (*fPcstream)<<"dcs"<< // Pulser information
929 "PulserTSector.="<<&sec<< // sector numbers
930 // // fit in respect to reference
931 "PulserTRef0.="<<&vec0<< // offset change
932 "PulserTRefY.="<<&vecLy<< // slope y change - rad
933 "PulserTRefX.="<<&vecLx<< // slope x change - rad
934 "PulserTRefChi2.="<<&vecChi2<< // chi2 (rms in cm)
935 "PulserTRefN.="<<&vecN<< //number of accepted points
936 // // fit in respect per mean per side
937 "PulserT0.="<<&vecA0<< // offset change
938 "PulserTY.="<<&vecALy<< // slope y change - rad
939 "PulserTX.="<<&vecALx<< // slope x change - rad
940 "PulserTChi2.="<<&vecAChi2; // chi2 (rms in cm)
941}
942
a3b590cf 943void AliTPCcalibSummary::ProcessCurrent(Int_t irun, Int_t itime){
944 //
945 // Dump current
946 //
947 //variables to export
948 //
949 static TObjArray *currentArray=new TObjArray(72); // current graphs
950 static TObjArray *currentArray2=new TObjArray(72); // current graphs to export
951 //
952 static TVectorD currentIROC(36); // current snapshots
953 static TVectorD currentOROC(36);
954 static TVectorF sector(72); //
955 static Double_t medcurIROC = 0;
956 static Double_t medcurOROC = 0;
957 //
958 static TVectorF minROC(72); // current mean +-5 minutes
959 static TVectorF maxROC(72);
960 static TVectorF meanROC(72);
961 static TVectorF medianROC(72);
962 static Double_t meanIIROC=0;
963 static Double_t meanIOROC=0;
964 static Double_t medianIIROC=0;
965 static Double_t medianIOROC=0;
966 //
967 AliDCSSensorArray* voltageArray = AliTPCcalibDB::Instance()->GetVoltageSensors(irun);
968 //
969 for(Int_t j=1; j<36; j++) currentIROC[j-1] = fCalibDB->GetChamberHighVoltage(irun, j,itime,-1,kTRUE);
970 for(Int_t j=36; j<72; j++) currentOROC[j-36] = fCalibDB->GetChamberHighVoltage(irun, j,itime,-1,kTRUE);
971 medcurIROC = TMath::Median(36, currentIROC.GetMatrixArray());
972 medcurOROC = TMath::Median(36, currentOROC.GetMatrixArray());
973
974
975 if (currentArray->At(0)==0){
976 for (Int_t isec=0; isec<72; isec++){
977 TString sensorName="";
978 const char* sideName=(isec%36<18) ? "A":"C";
979 if (isec<36){
980 //IROC
981 sensorName=Form("TPC_ANODE_I_%s%02d_IMEAS",sideName,isec%18);
982 }else{
983 //OROC
984 sensorName=Form("TPC_ANODE_O_%s%02d_0_IMEAS",sideName,isec%18);
985 }
986
987 AliDCSSensor *sensor = 0;
988 if (voltageArray) sensor= voltageArray->GetSensor(sensorName);
989 TGraph *gr=0;
990 if (!sensor) gr=new TGraph(1);
991 else{
992 if (!sensor->GetGraph()) gr=new TGraph(1);
993 else{
994 gr=sensor->GetGraph();
995 Double_t startTime=sensor->GetStartTime();
996 Double_t * time = new Double_t[gr->GetN()];
997 for (Int_t ip=0; ip<gr->GetN(); ip++){ time[ip]= (gr->GetX()[ip]*3600.)+startTime;}
998 gr=new TGraph(gr->GetN(), time, gr->GetY());
999 delete [] time;
1000 }
1001 }
1002 gr->Sort();
1003 currentArray->AddAt(gr, isec);
1004 currentArray->AddAt(gr->Clone(), isec);
1005 }
1006 }
1007
1008
1009 for (Int_t isec=0; isec<72; isec++){
1010 sector[isec]=isec;
1011 TGraph * gr = (TGraph*)currentArray->At(isec);
1012 TGraph * graph2 = (TGraph*)currentArray2->At(isec);
1013 Int_t firstBin= TMath::BinarySearch(gr->GetN(), gr->GetX(), itime-300.)-2;
1014 Int_t lastBin= TMath::BinarySearch(gr->GetN(), gr->GetX(), itime+300.)+2;
1015 if (firstBin<0) firstBin=0;
1016 if (lastBin>=gr->GetN()) lastBin=gr->GetN()-1;
1017 //
1018 if (firstBin<lastBin){
1019 //
1020 minROC[isec]=TMath::MinElement(lastBin-firstBin, &(gr->GetY()[firstBin]));
1021 maxROC[isec]=TMath::MaxElement(lastBin-firstBin, &(gr->GetY()[firstBin]));
1022 meanROC[isec]=TMath::Mean(lastBin-firstBin, &(gr->GetY()[firstBin]));
1023 medianROC[isec]=TMath::Median(lastBin-firstBin, &(gr->GetY()[firstBin]));
1024 graph2 = new TGraph(lastBin-firstBin, &(gr->GetX()[firstBin]), &(gr->GetY()[firstBin]));
1025 delete currentArray2->At(isec);
1026 currentArray2->AddAt(graph2,isec);
1027 }
1028 (*fPcstream)<<"dcs"<< // current information
1029 Form("current%d.=",isec)<<graph2;
1030 }
1031 meanIIROC=TMath::Mean(36, &(meanROC.GetMatrixArray()[0]));
1032 meanIOROC=TMath::Mean(36, &(meanROC.GetMatrixArray()[36]));
1033 medianIIROC=TMath::Median(36, &(meanROC.GetMatrixArray()[0]));
1034 medianIOROC=TMath::Median(36, &(meanROC.GetMatrixArray()[36]));
1035 //
1036 (*fPcstream)<<"dcs"<< // current information
1037 "isec.="<<&sector<< //sector number
1038 "IIROC.="<<&currentIROC<< // current sample at given moment
1039 "IOROC.="<<&currentOROC<< // current sample at given moment
1040 "medianIIROC="<<medcurIROC<< // median at given moment
1041 "medianIOROC="<<medcurOROC<< // median at given moment
1042 //
1043 "minIROC.="<<&minROC<< // minimum in +-5 min
1044 "maxIROC.="<<&maxROC<< // maximum in +-5 min
1045 "meanIROC.="<<&meanROC<< // mean in +-5 min
1046 "medianIROC.="<<&medianROC<< // median in +-5 min
1047 "meanIIROC5="<<meanIIROC<< // mean current in IROC +-5 minutes
1048 "meanIOROC5="<<meanIOROC<< // mean current in OROC
1049 "medianIIROC5="<<medianIIROC<< // median current in IROC
1050 "medianIOROC5="<<medianIOROC; // medianan current in OROC
1051
1052
1053 (*fPcstream)<<"current"<< // current information
1054 "time="<<itime<<
1055 "isec.="<<&sector<< //sector number
1056 "IIROC.="<<&currentIROC<< // current sample at given moment
1057 "IOROC.="<<&currentOROC<< // current sample at given moment
1058 "medianIIROC="<<medcurIROC<< // median at given moment
1059 "medianIOROC="<<medcurOROC<< // median at given moment
1060 //
1061 "minIROC.="<<&minROC<< // minimum in +-5 min
1062 "maxIROC.="<<&maxROC<< // maximum in +-5 min
1063 "meanIROC.="<<&meanROC<< // mean in +-5 min
1064 "medianIROC.="<<&medianROC<< // median in +-5 min
1065 "meanIIROC5="<<meanIIROC<< // mean current in IROC +-5 minutes
1066 "meanIOROC5="<<meanIOROC<< // mean current in OROC
1067 "medianIIROC5="<<medianIIROC<< // median current in IROC
1068 "medianIOROC5="<<medianIOROC<< // medianan current in OROC
1069 "\n";
1070
1071}
af6a50bb 1072
1073
1074
1075
1076// TCanvas * DrawCEDiff(TTree * tree){
1077
1078// TCanvas *canvasIO = new TCanvas("canvasCEIO","canvasCEIO");
1079// canvasIO->Divide(6,6);
1080// for (Int_t isec=0; isec<36; isec++){
1081// canvasIO->cd(isec+1);
1082// dcs->Draw(Form("CET0.fElements[%d]-CET0.fElements[%d]",isec+36,isec),Form("abs(CETRef0.fElements[%d])<0.3",isec),"");
1083// printf("%d\t%f\t%f\n",isec,dcs->GetHistogram()->GetMean(),dcs->GetHistogram()->GetRMS());
1084// }
1085
1086// }