]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/CalibMacros/CalibEnv.C
--This li
[u/mrichter/AliRoot.git] / TPC / CalibMacros / CalibEnv.C
CommitLineData
d47d42ee 1/*
7ce80437 2.x ~/NimStyle.C
3.x ~/rootlogon.C
4
d47d42ee 5gSystem->AddIncludePath("-I$ALICE_ROOT/TPC");
e2914767 6
e2914767 7
8
7ce80437 9.L $ALICE_ROOT/TPC/CalibMacros/CalibEnv.C+
d47d42ee 10Init();
e2914767 11CalibEnv("listAll.txt");
7ce80437 12GetTree();
13
e79211e8 14TFile f("dcsTime.root")
d47d42ee 15
7ce80437 16//
17// if you want to use alien OCDB
18//
19gSystem->Load("libXrdClient.so");
20gSystem->Load("libNetx.so");
21if (!gGrid) TGrid::Connect("alien://",0,0,"t");
22
23
e2914767 24
d47d42ee 25*/
26
27#include <iostream>
28#include <fstream>
5312f439 29#include <stdio.h>
d47d42ee 30#include <AliCDBManager.h>
31#include <AliCDBEntry.h>
32#include <AliLog.h>
33#include <AliMagF.h>
34#include "AliTPCcalibDB.h"
35#include "AliTPCAltroMapping.h"
36#include "AliTPCExB.h"
37#include "AliTPCCalROC.h"
38#include "AliTPCCalPad.h"
39#include "AliTPCSensorTempArray.h"
40#include "AliGRPObject.h"
41#include "AliTPCTransform.h"
42#include "TFile.h"
43#include "TKey.h"
44#include "TObjArray.h"
45#include "TObjString.h"
46#include "TString.h"
47#include "AliTPCCalPad.h"
5312f439 48#include "AliTPCROC.h"
49#include "AliTPCParam.h"
d47d42ee 50#include "AliTPCCalibPulser.h"
51#include "AliTPCCalibPedestal.h"
52#include "AliTPCCalibCE.h"
53#include "AliTPCExBFirst.h"
54#include "TTreeStream.h"
55#include "AliTPCTempMap.h"
5312f439 56#include "TVectorD.h"
57#include "TMatrixD.h"
d47d42ee 58
7ce80437 59
60TTree * dcsTree=0;
61
d47d42ee 62void ProcessGoofie( AliDCSSensorArray* goofieArray, TVectorD & vecEntries, TVectorD & vecMedian, TVectorD &vecMean, TVectorD &vecRMS);
5312f439 63void ProcessCEdata(const char* fitFormula, TVectorD &fitResultsA, TVectorD &fitResultsC);
64void ProcessNoiseData(TVectorD &vNoiseMean, TVectorD &vNoiseMeanSenRegions,
65 TVectorD &vNoiseRMS, TVectorD &vNoiseRMSSenRegions,
66 Int_t &nonMaskedZero);
67void ProcessPulser(Int_t &nMasked, Int_t &nonMaskedZero);
68void GetProductionInfo(Int_t run, Int_t &nalien, Int_t &nlocal);
69
d47d42ee 70void Init(){
71 //
72 //
73 //
162637e4 74 AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
75 AliCDBManager::Instance()->SetSpecificStorage("TPC/Calib/Parameters","local://$ALICE_ROOT/OCDB");
7ce80437 76 AliCDBManager::Instance()->SetSpecificStorage("GRP/GRP/Data","local:///lustre/alice/alien/alice/data/2008/LHC08d/OCDB/");
77 AliCDBManager::Instance()->SetSpecificStorage("TPC/Calib/Temperature","local:///lustre/alice/alien/alice/data/2008/LHC08d/OCDB/");
78 AliCDBManager::Instance()->SetSpecificStorage("TPC/Calib/HighVoltage","local:///lustre/alice/alien/alice/data/2008/LHC08d/OCDB/");
79 AliCDBManager::Instance()->SetSpecificStorage("TPC/Calib/Goofie","local:///lustre/alice/alien/alice/data/2008/LHC08d/OCDB/");
d47d42ee 80 AliCDBManager::Instance()->SetRun(1);
e2914767 81}
82
d47d42ee 83
e2914767 84void InitAlien(const char *path="LHC08b"){
85 //
86 //
87 //
88 TString alpath="alien://folder=/alice/data/2008/";
89 alpath+=path;
90 alpath+="/OCDB";
7ce80437 91
162637e4 92 AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
93 AliCDBManager::Instance()->SetSpecificStorage("TPC/Calib/Parameters","local://$ALICE_ROOT/OCDB");
e2914767 94 AliCDBManager::Instance()->SetSpecificStorage("GRP/GRP/Data",alpath.Data());
95 AliCDBManager::Instance()->SetSpecificStorage("TPC/Calib/Temperature",alpath.Data());
96 AliCDBManager::Instance()->SetSpecificStorage("TPC/Calib/Goofie",alpath.Data());
97 AliCDBManager::Instance()->SetRun(1);
d47d42ee 98}
99
100
101void CalibEnv(const char * runList){
102 //
103 //
7ce80437 104 //
e79211e8 105 AliTPCcalibDB * calibDB = AliTPCcalibDB::Instance();
d47d42ee 106 ifstream in;
107 in.open(runList);
108 Int_t irun=0;
e79211e8 109 TTreeSRedirector *pcstream = new TTreeSRedirector("dcsTime.root");
d47d42ee 110 // for (Int_t irun=startRun; irun<stopRun; irun++){
111 while(in.good()) {
112 in >> irun;
113 if (irun==0) continue;
e2914767 114 printf("Processing run %d ...\n",irun);
5312f439 115 AliTPCcalibDB::Instance()->SetRun(irun);
d47d42ee 116 AliDCSSensor * sensorPressure = AliTPCcalibDB::Instance()->GetPressureSensor(irun);
117 if (!sensorPressure) continue;
118 AliTPCSensorTempArray * tempArray = AliTPCcalibDB::Instance()->GetTemperatureSensor(irun);
119 AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
120 AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(irun);
121 //
122 Int_t startTime = sensorPressure->GetStartTime();
123 Int_t endTime = sensorPressure->GetEndTime();
124 Int_t dtime = TMath::Max((endTime-startTime)/20,10*60);
5312f439 125 //CE data processing - see ProcessCEdata function for description of the results
126 TVectorD fitResultsA, fitResultsC;
127 ProcessCEdata("gx++gy++lx++lx**2",fitResultsA,fitResultsC);
128 //noise data Processing - see ProcessNoiseData function for description of the results
129 TVectorD vNoiseMean, vNoiseMeanSenRegions, vNoiseRMS, vNoiseRMSSenRegions;
130 Int_t nonMaskedZero=0;
131 ProcessNoiseData(vNoiseMean, vNoiseMeanSenRegions, vNoiseRMS, vNoiseRMSSenRegions, nonMaskedZero);
132 //L3 data
133 Float_t bz=AliTPCcalibDB::GetBz(irun);
134 Char_t l3pol=AliTPCcalibDB::GetL3Polarity(irun);
135 //calibration Pulser data processing
136 Int_t nMasked=0;
137 Int_t nOffChannels=0;
138 ProcessPulser(nMasked,nOffChannels);
139 //production information
140 Int_t nalien=0,nlocal=0;
141 GetProductionInfo(irun, nalien, nlocal);
142 //run type
143 TObjString runType(AliTPCcalibDB::GetRunType(irun).Data());
144
d47d42ee 145 for (Int_t itime=startTime; itime<endTime; itime+=dtime){
146 //
5312f439 147 TTimeStamp tstamp(itime);
e79211e8 148 Float_t valuePressure = calibDB->GetPressure(tstamp,irun,0);
149 Float_t valuePressure2 = calibDB->GetPressure(tstamp,irun,1);
5312f439 150 //temperature fits
d47d42ee 151 TLinearFitter * fitter = 0;
152 TVectorD vecTemp[10];
e79211e8 153 for (Int_t itype=0; itype<5; itype++)
5312f439 154 for (Int_t iside=0; iside<2; iside++){
155 fitter= tempMap->GetLinearFitter(itype,iside,tstamp);
156 if (!fitter) continue;
157 fitter->Eval();
158 fitter->GetParameters(vecTemp[itype+iside*5]);
159 delete fitter;
160 }
161 //measured skirt temperatures
162 TVectorD vecSkirtTempA(18);
163 TVectorD vecSkirtTempC(18);
164 Int_t nsenTemp=tempArray->NumSensors();
165 for (Int_t isenTemp=0;isenTemp<nsenTemp;++isenTemp){
166 AliTPCSensorTemp *senTemp=(AliTPCSensorTemp*)tempArray->GetSensorNum(isenTemp);
167 if (senTemp->GetType()!=3) continue;
168 if (TMath::Sqrt(senTemp->GetX()*senTemp->GetX()+senTemp->GetY()*senTemp->GetY())<100) continue; //only skirt, outer FC vessel
169 Double_t val=senTemp->GetValue(tstamp);
170 if (senTemp->GetSide()==0)
171 vecSkirtTempA[senTemp->GetSector()]=val;
172 else
173 vecSkirtTempC[senTemp->GetSector()]=val;
174 }
d47d42ee 175
176 TVectorD vecGoofie, vecEntries, vecMean, vecMedian,vecRMS;
5312f439 177 if (goofieArray){
178 vecGoofie.ResizeTo(goofieArray->NumSensors());
179 ProcessGoofie(goofieArray, vecEntries ,vecMedian, vecMean, vecRMS);
180 //
181 for (Int_t isensor=0; isensor<goofieArray->NumSensors();isensor++){
182 AliDCSSensor *gsensor = goofieArray->GetSensor(isensor);
183 if (gsensor){
184 vecGoofie[isensor] = gsensor->GetValue(tstamp);
185 }
186 }
d47d42ee 187 }
e79211e8 188 Double_t ptrelative0 = AliTPCcalibDB::GetPTRelative(tstamp,irun,0);
189 Double_t ptrelative1 = AliTPCcalibDB::GetPTRelative(tstamp,irun,1);
e2914767 190 //
5312f439 191 TVectorD voltagesIROC(36);
192 TVectorD voltagesOROC(36);
7ce80437 193 for(Int_t j=1; j<36; j++) voltagesIROC[j-1] = AliTPCcalibDB::Instance()->GetChamberHighVoltage(irun, j,tstamp);
194 for(Int_t j=36; j<72; j++) voltagesOROC[j-36] = AliTPCcalibDB::Instance()->GetChamberHighVoltage(irun, j,tstamp);
195 Double_t voltIROC = TMath::Median(36, voltagesIROC.GetMatrixArray());
196 Double_t voltOROC = TMath::Median(36, voltagesOROC.GetMatrixArray());
197 //
5312f439 198 Float_t coverIA=AliTPCcalibDB::GetCoverVoltage(irun,0,itime);
199 Float_t coverIC=AliTPCcalibDB::GetCoverVoltage(irun,18,itime);
200 Float_t coverOA=AliTPCcalibDB::GetCoverVoltage(irun,36,itime);
201 Float_t coverOC=AliTPCcalibDB::GetCoverVoltage(irun,54,itime);
202 Float_t skirtA=AliTPCcalibDB::GetSkirtVoltage(irun,0,itime);
203 Float_t skirtC=AliTPCcalibDB::GetSkirtVoltage(irun,18,itime);
204 Float_t ggOffA=AliTPCcalibDB::GetGGoffsetVoltage(irun,0,itime);
205 Float_t ggOffC=AliTPCcalibDB::GetGGoffsetVoltage(irun,18,itime);
206
207
208
e2914767 209 //tempMap->GetLinearFitter(0,0,itime);
d47d42ee 210 (*pcstream)<<"dcs"<<
5312f439 211 "run="<<irun<<
212 "time="<<itime<<
213 //run type
214 "runType.="<<&runType<<
215 // voltage setting
216 "VIROC.="<<&voltagesIROC<<
217 "VOROC.="<<&voltagesOROC<<
218 "medianVIROC="<<voltIROC<<
219 "medianVOROC="<<voltOROC<<
7ce80437 220 "coverIA=" << coverIA <<
221 "coverIC=" << coverIC <<
222 "coverOA=" << coverOA <<
223 "coverOC=" << coverOC <<
224 "skirtA=" << skirtA <<
225 "skirtC=" << skirtC <<
226 "ggOffA=" << ggOffA <<
227 "ggOffC=" << ggOffC <<
5312f439 228 //
229 "ptrel0="<<ptrelative0<< // deltaTP/TP - A side
230 "ptrel1="<<ptrelative1<< // deltaTP/TPC - C side
231 "goofie.="<<&vecGoofie<<
232 "goofieE.="<<&vecEntries<<
233 "goofieMean.="<<&vecMean<<
234 "goofieMedian.="<<&vecMedian<<
235 "goofieRMS.="<<&vecRMS<<
236 //
237 "press="<<valuePressure<<
238 "press2="<<valuePressure2<<
239 "temp00.="<<&vecTemp[0]<<
240 "temp10.="<<&vecTemp[1]<<
241 "temp20.="<<&vecTemp[2]<<
242 "temp30.="<<&vecTemp[3]<<
243 "temp40.="<<&vecTemp[4]<<
244 "temp01.="<<&vecTemp[5]<<
245 "temp11.="<<&vecTemp[6]<<
246 "temp21.="<<&vecTemp[7]<<
247 "temp31.="<<&vecTemp[8]<<
248 "temp41.="<<&vecTemp[9]<<
249 "tempSkirtA.="<<&vecSkirtTempA<<
250 "tempSkirtC.="<<&vecSkirtTempC<<
251 //noise data
252 "meanNoise.="<<&vNoiseMean<<
253 "meanNoiseSen.="<<&vNoiseMeanSenRegions<<
254 "rmsNoise.="<<&vNoiseRMS<<
255 "rmsNoiseSen.="<<&vNoiseRMSSenRegions<<
256 "zeroNoise="<<nonMaskedZero<<
257 //pulser data
258 "nMasked="<<nMasked<< //should perhaps go to altro data
259 "nOffPulser="<<nOffChannels<<
260 //ce data
261 "CEfitA.="<<&fitResultsA<<
262 "CEfitC.="<<&fitResultsC<<
263 // b field
264 "Bz="<< bz <<
265 "L3polarity="<<l3pol<<
266 // production information
267 "nalien="<<nalien<<
268 "nlocal="<<nlocal<<
269 "\n";
d47d42ee 270 }
271 }
272 delete pcstream;
273}
274
275
276void ProcessGoofie( AliDCSSensorArray* goofieArray, TVectorD & vecEntries, TVectorD & vecMedian, TVectorD &vecMean, TVectorD &vecRMS){
277 /*
7ce80437 278
d47d42ee 279 1 TPC_ANODE_I_A00_STAT
280 2 TPC_DVM_CO2
281 3 TPC_DVM_DriftVelocity
282 4 TPC_DVM_FCageHV
283 5 TPC_DVM_GainFar
284 6 TPC_DVM_GainNear
285 7 TPC_DVM_N2
286 8 TPC_DVM_NumberOfSparks
287 9 TPC_DVM_PeakAreaFar
288 10 TPC_DVM_PeakAreaNear
289 11 TPC_DVM_PeakPosFar
290 12 TPC_DVM_PeakPosNear
291 13 TPC_DVM_PickupHV
292 14 TPC_DVM_Pressure
293 15 TPC_DVM_T1_Over_P
294 16 TPC_DVM_T2_Over_P
295 17 TPC_DVM_T_Over_P
296 18 TPC_DVM_TemperatureS1
297 */
298 //
299 //
300 // TVectorD vecMedian; TVectorD vecEntries; TVectorD vecMean; TVectorD vecRMS;
301 Double_t kEpsilon=0.0000000001;
302 Double_t kBig=100000000000.;
303 Int_t nsensors = goofieArray->NumSensors();
304 vecEntries.ResizeTo(nsensors);
305 vecMedian.ResizeTo(nsensors);
306 vecMean.ResizeTo(nsensors);
307 vecRMS.ResizeTo(nsensors);
308 TVectorF values;
309 for (Int_t isensor=0; isensor<goofieArray->NumSensors();isensor++){
310 AliDCSSensor *gsensor = goofieArray->GetSensor(isensor);
311 if (gsensor && gsensor->GetGraph()){
312 Int_t npoints = gsensor->GetGraph()->GetN();
313 // filter zeroes
314 values.ResizeTo(npoints);
315 Int_t nused =0;
316 for (Int_t ipoint=0; ipoint<npoints; ipoint++){
7ce80437 317 if (TMath::Abs(gsensor->GetGraph()->GetY()[ipoint])>kEpsilon &&
318 TMath::Abs(gsensor->GetGraph()->GetY()[ipoint])<kBig ){
319 values[nused]=gsensor->GetGraph()->GetY()[ipoint];
320 nused++;
321 }
d47d42ee 322 }
323 //
7ce80437 324 vecEntries[isensor]= nused;
d47d42ee 325 if (nused>1){
7ce80437 326 vecMedian[isensor] = TMath::Median(nused,values.GetMatrixArray());
327 vecMean[isensor] = TMath::Mean(nused,values.GetMatrixArray());
328 vecRMS[isensor] = TMath::RMS(nused,values.GetMatrixArray());
d47d42ee 329 }
330 }
331 }
332}
333
5312f439 334void ProcessCEdata(const char* fitFormula, TVectorD &fitResultsA, TVectorD &fitResultsC)
335{
336 //
337 // Process the CE data for this run
338 // the return TVectorD arrays contian the results of the fit
339 //
340 const Float_t irocToffset=0.2;
341 const Float_t tMaxLimit=1.2;
342 //retrieve CE and ALTRO data
343 AliTPCCalPad *cet0=AliTPCcalibDB::Instance()->GetCETmean();
344 if (!cet0){
345 TString fitString(fitFormula);
346 fitString.ReplaceAll("++","#");
347 Int_t ndim=fitString.CountChar('#')+1;
348 fitResultsA.ResizeTo(ndim);
349 fitResultsC.ResizeTo(ndim);
350 return;
351 }
352 AliTPCCalPad padT0(*cet0);
353 AliTPCCalPad *padMasked=AliTPCcalibDB::Instance()->GetALTROMasked();
354 //create outlier map
355 AliTPCCalPad out("out","out");
356 //loop over all channels
357 for (UInt_t iroc=0;iroc<padT0.kNsec;++iroc){
358 AliTPCCalROC *rocData=padT0.GetCalROC(iroc);
359 AliTPCCalROC *rocMasked=padMasked->GetCalROC(iroc);
360 AliTPCCalROC *rocOut=out.GetCalROC(iroc);
361 if (!rocData) continue;
362 //add time offset to IROCs
363 if (iroc<AliTPCROC::Instance()->GetNInnerSector())
364 rocData->Add(irocToffset);
365 //select outliers
366 for (UInt_t ichannel=0;ichannel<rocData->GetNchannels();++ichannel){
367 if (rocMasked && rocMasked->GetValue(ichannel)) rocOut->SetValue(ichannel,1);
368 Float_t valTmean=rocData->GetValue(ichannel);
369 if (valTmean==0) rocOut->SetValue(ichannel,1); //exclude values that are exactly 0
370 if (TMath::Abs(valTmean)>tMaxLimit) rocOut->SetValue(ichannel,1); // exclude channels with too large variations
371 }
372 }
373 //perform fit
374 TMatrixD dummy;
375 Float_t chi2A,chi2C;
376 padT0.GlobalSidesFit(&out,fitFormula,fitResultsA,fitResultsC,dummy,dummy,chi2A,chi2C);
377}
378
379void ProcessNoiseData(TVectorD &vNoiseMean, TVectorD &vNoiseMeanSenRegions,
380 TVectorD &vNoiseRMS, TVectorD &vNoiseRMSSenRegions,
381 Int_t &nonMaskedZero)
382{
383 //
384 // process noise data
385 // vNoiseMean/RMS contains the Mean/RMS noise of the complete TPC [0], IROCs only [1],
386 // OROCs small pads [2] and OROCs large pads [3]
387 // vNoiseMean/RMSsenRegions constains the same information, but only for the sensitive regions (edge pads, corners, IROC spot)
388 //
389
390 //set proper size and reset
391 const UInt_t infoSize=4;
392 vNoiseMean.ResizeTo(infoSize);
393 vNoiseMeanSenRegions.ResizeTo(infoSize);
394 vNoiseRMS.ResizeTo(infoSize);
395 vNoiseRMSSenRegions.ResizeTo(infoSize);
396 vNoiseMean.Zero();
397 vNoiseMeanSenRegions.Zero();
398 vNoiseRMS.Zero();
399 vNoiseRMSSenRegions.Zero();
400 //counters
401 TVectorD c(infoSize);
402 TVectorD cs(infoSize);
403 //tpc parameters
404 AliTPCParam par;
405 par.Update();
406 //retrieve noise and ALTRO data
407 AliTPCCalPad *noise=AliTPCcalibDB::Instance()->GetPadNoise();
408 AliTPCCalPad *padMasked=AliTPCcalibDB::Instance()->GetALTROMasked();
409 //create IROC, OROC1, OROC2 and sensitive region masks
410 for (UInt_t isec=0;isec<AliTPCCalPad::kNsec;++isec){
411 AliTPCCalROC *noiseROC=noise->GetCalROC(isec);
412 UInt_t nrows=noiseROC->GetNrows();
413 for (UInt_t irow=0;irow<nrows;++irow){
414 UInt_t npads=noiseROC->GetNPads(irow);
415 for (UInt_t ipad=0;ipad<npads;++ipad){
416 Int_t masked=(Int_t)padMasked->GetCalROC(isec)->GetValue(irow,ipad);
417 Float_t noiseVal=noiseROC->GetValue(irow,ipad);
418 if (masked) continue; // don't use inactive pads
419 //check if noise==0
420 if (noiseVal==0) ++nonMaskedZero;
421 Int_t cpad=(Int_t)ipad-(Int_t)npads/2;
422 Int_t masksen=1; // sensitive pards are not masked (0)
423 if (ipad<2||npads-ipad-1<2) masksen=0; //don't mask edge pads (sensitive)
424 if (isec<AliTPCROC::Instance()->GetNInnerSector()){
425 //IROCs
426 if (irow>19&&irow<46){
427 if (TMath::Abs(cpad)<7) masksen=0; //IROC spot
428 }
429 Int_t type=1;
430 vNoiseMean[type]+=noiseVal;
431 vNoiseRMS[type]+=noiseVal*noiseVal;
432 ++c[type];
433 if (!masksen){
434 vNoiseMeanSenRegions[type]+=noiseVal;
435 vNoiseRMSSenRegions[type]+=noiseVal*noiseVal;
436 ++cs[type];
437 }
438 } else {
439 //OROCs
440 //define sensive regions
441 if ((nrows-irow-1)<3) masksen=0; //last three rows in OROCs are sensitive
442 if ( irow>75 ){
443 Int_t padEdge=(Int_t)TMath::Min(ipad,npads-ipad);
444 if (padEdge<((((Int_t)irow-76)/4+1))*2) masksen=0; //OROC outer corners are sensitive
445 }
446 if ((Int_t)irow<par.GetNRowUp1()){
447 //OROC1
448 Int_t type=2;
449 vNoiseMean[type]+=noiseVal;
450 vNoiseRMS[type]+=noiseVal*noiseVal;
451 ++c[type];
452 if (!masksen){
453 vNoiseMeanSenRegions[type]+=noiseVal;
454 vNoiseRMSSenRegions[type]+=noiseVal*noiseVal;
455 ++cs[type];
456 }
457 }else{
458 //OROC2
459 Int_t type=3;
460 vNoiseMean[type]+=noiseVal;
461 vNoiseRMS[type]+=noiseVal*noiseVal;
462 ++c[type];
463 if (!masksen){
464 vNoiseMeanSenRegions[type]+=noiseVal;
465 vNoiseRMSSenRegions[type]+=noiseVal*noiseVal;
466 ++cs[type];
467 }
468 }
469 }
470 //whole tpc
471 Int_t type=0;
472 vNoiseMean[type]+=noiseVal;
473 vNoiseRMS[type]+=noiseVal*noiseVal;
474 ++c[type];
475 if (!masksen){
476 vNoiseMeanSenRegions[type]+=noiseVal;
477 vNoiseRMSSenRegions[type]+=noiseVal*noiseVal;
478 ++cs[type];
479 }
480 }//end loop pads
481 }//end loop rows
482 }//end loop sectors (rocs)
483
484 //calculate mean and RMS
485 const Double_t verySmall=0.0000000001;
486 for (UInt_t i=0;i<infoSize;++i){
487 Double_t mean=0;
488 Double_t rms=0;
489 Double_t meanSen=0;
490 Double_t rmsSen=0;
491
492 if (c[i]>verySmall){
493 mean=vNoiseMean[i]/c[i];
494 rms=vNoiseRMS[i];
495 rms=TMath::Sqrt(TMath::Abs(rms/c[i]-mean*mean));
496 }
497 vNoiseMean[i]=mean;
498 vNoiseRMS[i]=rms;
499
500 if (cs[i]>verySmall){
501 meanSen=vNoiseMeanSenRegions[i]/cs[i];
502 rmsSen=vNoiseRMSSenRegions[i];
503 rmsSen=TMath::Sqrt(TMath::Abs(rmsSen/cs[i]-meanSen*meanSen));
504 }
505 vNoiseMeanSenRegions[i]=meanSen;
506 vNoiseRMSSenRegions[i]=rmsSen;
507 }
508}
509
510void ProcessPulser(Int_t &nMasked, Int_t &nonMaskedZero)
511{
512 //
513 // Process the Pulser information
514 //
515
516 //reset counters
517 nonMaskedZero=0;
518 nMasked=0;
519 //retrieve pulser and ALTRO data
520 AliTPCCalPad *pulserTmean=AliTPCcalibDB::Instance()->GetPulserTmean();
521// AliTPCCalPad *pulserTrms=AliTPCcalibDB::Instance()->GetPulserTrms();
522// AliTPCCalPad *pulserQmean=AliTPCcalibDB::Instance()->GetPulserQmean();
523 AliTPCCalPad *padMasked=AliTPCcalibDB::Instance()->GetALTROMasked();
524 //create IROC, OROC1, OROC2 and sensitive region masks
525 for (UInt_t isec=0;isec<AliTPCCalPad::kNsec;++isec){
526 AliTPCCalROC *tmeanROC=pulserTmean->GetCalROC(isec);
527// AliTPCCalROC *trmsROC=pulserTrms->GetCalROC(isec);
528// AliTPCCalROC *qmeanROC=pulserQmean->GetCalROC(isec);
529 Float_t tmeanMedian=tmeanROC->GetMedian();
530 UInt_t nrows=tmeanROC->GetNrows();
531 for (UInt_t irow=0;irow<nrows;++irow){
532 UInt_t npads=tmeanROC->GetNPads(irow);
533 for (UInt_t ipad=0;ipad<npads;++ipad){
534 Int_t masked=(Int_t)padMasked->GetCalROC(isec)->GetValue(irow,ipad);
535 Float_t tmeanVal=tmeanROC->GetValue(irow,ipad);
536// Float_t trmsVal =trmsROC->GetValue(irow,ipad);
537// Float_t qmeanVal=qmeanROC->GetValue(irow,ipad);
538 if (masked){
539 ++nMasked;
540 continue; // don't use inactive pads
541 }
542 if ( TMath::Abs(tmeanVal-tmeanMedian)>1.5 ) ++nonMaskedZero;
543 }
544 }
545 }
546}
547
548void GetProductionInfo(Int_t run, Int_t &nalien, Int_t &nlocal){
549 //
550 // find number of ESDs in central and local production for run
551 //
552
553 nalien=0;
554 nlocal=0;
555 TString sNlines;
556 FILE *pipe = 0x0;
557 //find number of ESDs in alien
558 TString command="alien_find /alice/data/2009 ";
559 command += Form("%09d",run);
560 command += " | grep AliESDs.root | wc -l";
561 pipe=gSystem->OpenPipe(command.Data(),"r");
562 sNlines.Gets(pipe);
563 gSystem->ClosePipe(pipe);
564 nalien=sNlines.Atoi();
565 //find number of ESDs local
566 command="find /lustre/alice/alien/alice/data/2009 -name ";
567 command += Form("%09d",run);
568 command += " | grep AliESDs.root | wc -l";
569 pipe = gSystem->OpenPipe(command.Data(),"r");
570 sNlines.Gets(pipe);
571 gSystem->ClosePipe(pipe);
572 nlocal=sNlines.Atoi();
573}
e2914767 574
575void FilterMag(const char * runList){
576 //
577 //
578 //
579 // AliTPCcalibDB * calibDB = AliTPCcalibDB::Instance();
580 ifstream in;
581 in.open(runList);
582 Int_t irun=0;
583 while(in.good()) {
584 in >> irun;
585 if (irun==0) continue;
586 AliGRPObject *grp = AliTPCcalibDB::GetGRP(irun);
587 Float_t current = -1;
588 Float_t bz = -1;
5312f439 589// Float_t press = 0;
e2914767 590 if (grp){
591 current = grp->GetL3Current((AliGRPObject::Stats)0);
592 bz = 5*current/30000.;
593 printf("Run%d\tL3 current%f\tBz\t%f\n",irun,current,bz);
594 }
595 else{
596 printf("Run%d\tL3 current%f\tBz\t%f\n",irun,current,bz);
597 }
598 }
7ce80437 599
600}
601
602
603void GetTree(){
604 TFile *fdcs = new TFile("dcsTime.root");
5312f439 605 dcsTree = (TTree*)fdcs->Get("dcs");
7ce80437 606 //
607 // mean temp A
5312f439 608
7ce80437 609 dcsTree->Draw("temp30.fElements[0]");
5312f439 610
e2914767 611}
612
7ce80437 613void GetNominalValues(){
614 //
615 if (!dcsTree) return;
616}
617
618
619
620
d47d42ee 621/*
622
623AliDCSSensor * sensorPressure = AliTPCcalibDB::Instance()->GetPressureSensor(62084);
624entry = AliCDBManager::Instance()->Get("TPC/Calib/Temperature",run);
625AliTPCSensorTempArray * tempArray = (AliTPCSensorTempArray *)entry->GetObject();
626AliTPCSensorTempArray * tempArray = (AliTPCSensorTempArray *)AliTPCcalibDB::Instance()->GetTemperatureSensor(62084)
627AliTPCTempMap * tempMap = new AliTPCTempMap(tempArray);
e79211e8 628TLinearFitter * fitter = tempMap->GetLinearFitter(0,0,tempArray->GetStartTime());
d47d42ee 629
630AliDCSSensorArray* goofieArray = AliTPCcalibDB::Instance()->GetGoofieSensors(62084);
631
632*/
7ce80437 633
634/*
635
636void PlotPressureResol(){
637 //
638 // Example
639 //
640 dcs->Draw("100000*(press-press2-4.782)/press/sqrt(2.)>>his(100,-50,50)","run>61400","")
641 his->SetXTitle("#sigma_{P/P_{0}}(10^{-5})");
642 gPad->SaveAs("picDCS/deltaPoverP.eps");
643 gPad->SaveAs("picDCS/deltaPoverP.gif");
5312f439 644
7ce80437 645}
646void PlotTresol(){
647 //
648 // T resolution example
649 // plot difference of the temperature from A and C side
650 // Supposing the error is independent - (division by sqrt(2))
651 dcs->Draw("100000*(temp30.fElements[0]-temp31.fElements[0]+0.00509)/(temp31.fElements[0]+273.15)/sqrt(2.)>>his(100,-5,5)","run>61400","");
652 his->SetXTitle("#sigma_{T/T_{0}}(10^{-5})");
653 gPad->SaveAs("picDCS/deltaToverT.eps");
654 gPad->SaveAs("picDCS/deltaToverT.gif");
655}
656*/