1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
19 Responsible: Raphaelle Bailhache (rbailhache@ikf.uni-frankfurt.de)
20 Code to analyze the TRD calibration and to produce OCDB entries
24 gSystem->Load("libANALYSIS");
25 gSystem->Load("libTRDcalib");
27 AliTRDPreprocessorOffline proces;
28 TString ocdbPath="local:////"
29 ocdbPath+=gSystem->GetFromPipe("pwd");
31 proces.CalibTimeGain("CalibObjects.root",run0,run1,ocdbPath);
32 proces.CalibTimeVdrift("CalibObjects.root",run0,run1,ocdbPath);
33 // take the raw calibration data from the file CalibObjects.root
34 // and make a OCDB entry with run validity run0-run1
35 // results are stored at the ocdbPath - local or alien ...
36 // default storage ""- data stored at current working directory
40 #include "Riostream.h"
49 #include "TProfile2D.h"
50 #include "AliTRDCalDet.h"
51 #include "AliTRDCalPad.h"
52 #include "AliCDBMetaData.h"
54 #include "AliCDBManager.h"
55 #include "AliCDBStorage.h"
56 #include "AliTRDCalibraMode.h"
57 #include "AliTRDCalibraFit.h"
58 #include "AliTRDCalibraVdriftLinearFit.h"
59 #include "AliTRDPreprocessorOffline.h"
60 #include "AliTRDCalChamberStatus.h"
63 ClassImp(AliTRDPreprocessorOffline)
65 AliTRDPreprocessorOffline::AliTRDPreprocessorOffline():
66 TNamed("TPCPreprocessorOffline","TPCPreprocessorOffline"),
68 fNameList("TRDCalib"),
70 fCalDetVdriftUsed(0x0),
74 fAliTRDCalibraVdriftLinearFit(0x0),
77 fPlots(new TObjArray(9)),
78 fCalibObjects(new TObjArray(9)),
80 fSubVersionGainUsed(0),
81 fFirstRunVdriftUsed(0),
82 fVersionVdriftUsed(0),
83 fSubVersionVdriftUsed(0),
84 fSwitchOnValidation(kTRUE),
85 fVdriftValidated(kFALSE),
87 fMinStatsVdriftT0PH(800*20),
88 fMinStatsVdriftLinear(800),
91 fBackCorrectGain(kFALSE),
92 fBackCorrectVdrift(kTRUE),
93 fNotEnoughStatisticsForTheGain(kFALSE),
94 fNotEnoughStatisticsForTheVdriftLinear(kFALSE),
98 // default constructor
101 //_________________________________________________________________________________________________________________
102 AliTRDPreprocessorOffline::~AliTRDPreprocessorOffline() {
107 if(fCalDetGainUsed) delete fCalDetGainUsed;
108 if(fCalDetVdriftUsed) delete fCalDetVdriftUsed;
109 if(fCH2d) delete fCH2d;
110 if(fPH2d) delete fPH2d;
111 if(fPRF2d) delete fPRF2d;
112 if(fAliTRDCalibraVdriftLinearFit) delete fAliTRDCalibraVdriftLinearFit;
113 if(fNEvents) delete fNEvents;
114 if(fAbsoluteGain) delete fAbsoluteGain;
115 if(fPlots) delete fPlots;
116 if(fCalibObjects) delete fCalibObjects;
119 //___________________________________________________________________________________________________________________
121 void AliTRDPreprocessorOffline::CalibVdriftT0(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){
123 // make calibration of the drift velocity
125 // file - the location of input file
126 // startRunNumber, endRunNumber - run validity period
127 // ocdbStorage - path to the OCDB storage
128 // - if empty - local storage 'pwd' uesed
129 if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
133 fVdriftValidated = kTRUE;
134 fT0Validated = kTRUE;
135 fNotEnoughStatisticsForTheVdriftLinear = kFALSE;
137 // 2. extraction of the information
139 if(ReadVdriftLinearFitGlobal(file)) AnalyzeVdriftLinearFit();
140 if(ReadVdriftT0Global(file)) AnalyzeVdriftT0();
142 // 3. Append QA plots
144 //MakeDefaultPlots(fVdriftArray,fVdriftArray);
147 // 4. validate OCDB entries
149 if(fSwitchOnValidation==kTRUE && ValidateVdrift()==kFALSE) {
150 AliError("TRD vdrift OCDB parameters out of range!");
151 fVdriftValidated = kFALSE;
153 if(fSwitchOnValidation==kTRUE && ValidateT0()==kFALSE) {
154 AliError("TRD t0 OCDB parameters out of range!");
155 fT0Validated = kFALSE;
161 if(fVdriftValidated) UpdateOCDBVdrift(startRunNumber,endRunNumber,ocdbStorage);
162 if(fT0Validated) UpdateOCDBT0(startRunNumber,endRunNumber,ocdbStorage);
163 UpdateOCDBExB(startRunNumber,endRunNumber,ocdbStorage);
166 //_________________________________________________________________________________________________________________
168 void AliTRDPreprocessorOffline::CalibGain(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){
170 // make calibration of the drift velocity
172 // file - the location of input file
173 // startRunNumber, endRunNumber - run validity period
174 // ocdbStorage - path to the OCDB storage
175 // - if empty - local storage 'pwd' uesed
176 if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
178 fNotEnoughStatisticsForTheGain = kFALSE;
181 if(!ReadGainGlobal(file)) return;
184 // 2. extraction of the information
187 if(fBackCorrectGain) CorrectFromDetGainUsed();
188 if(fBackCorrectVdrift) CorrectFromDetVdriftUsed();
190 // 3. Append QA plots
192 //MakeDefaultPlots(fVdriftArray,fVdriftArray);
195 // 4. validate OCDB entries
197 if(fSwitchOnValidation==kTRUE && ValidateGain()==kFALSE) {
198 AliError("TRD gain OCDB parameters out of range!");
205 if((!fCalDetVdriftUsed) || (fCalDetVdriftUsed && fVdriftValidated)) UpdateOCDBGain(startRunNumber,endRunNumber,ocdbStorage);
209 //________________________________________________________________________________________________________________
211 void AliTRDPreprocessorOffline::CalibPRF(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){
213 // make calibration of the drift velocity
215 // file - the location of input file
216 // startRunNumber, endRunNumber - run validity period
217 // ocdbStorage - path to the OCDB storage
218 // - if empty - local storage 'pwd' uesed
219 if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
222 if(!ReadPRFGlobal(file)) return;
225 // 2. extraction of the information
229 // 3. Append QA plots
231 //MakeDefaultPlots(fVdriftArray,fVdriftArray);
235 // 4. validate OCDB entries
237 if(fSwitchOnValidation==kTRUE && ValidatePRF()==kFALSE) {
238 AliError("TRD prf OCDB parameters out of range!");
245 UpdateOCDBPRF(startRunNumber,endRunNumber,ocdbStorage);
248 //________________________________________________________________________________________________________________
250 void AliTRDPreprocessorOffline::CalibChamberStatus(Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){
252 // make calibration of the chamber status
254 // startRunNumber, endRunNumber - run validity period
255 // ocdbStorage - path to the OCDB storage
256 // - if empty - local storage 'pwd' uesed
257 if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
260 // 2. extraction of the information
262 if(!AnalyzeChamberStatus()) return;
264 // 3. Append QA plots
266 //MakeDefaultPlots(fVdriftArray,fVdriftArray);
270 // 4. validate OCDB entries
272 if(fSwitchOnValidation==kTRUE && ValidateChamberStatus()==kFALSE) {
273 AliError("TRD Chamber status OCDB parameters not ok!");
280 if((!fNotEnoughStatisticsForTheGain) && (!fNotEnoughStatisticsForTheGain)) UpdateOCDBChamberStatus(startRunNumber,endRunNumber,ocdbStorage);
283 //______________________________________________________________________________________________________
284 Bool_t AliTRDPreprocessorOffline::Init(const Char_t* fileName){
286 // read the calibration used during the reconstruction
289 if(ReadVdriftT0Global(fileName)) {
291 TString nameph = fPH2d->GetTitle();
292 fFirstRunVdriftUsed = GetFirstRun(nameph);
293 fVersionVdriftUsed = GetVersion(nameph);
294 fSubVersionVdriftUsed = GetSubVersion(nameph);
296 //printf("Found Version %d, Subversion %d for vdrift\n",fVersionVdriftUsed,fSubVersionVdriftUsed);
300 if(ReadGainGlobal(fileName)) {
302 TString namech = fCH2d->GetTitle();
303 fVersionGainUsed = GetVersion(namech);
304 fSubVersionGainUsed = GetSubVersion(namech);
306 //printf("Found Version %d, Subversion %d for gain\n",fVersionGainUsed,fSubVersionGainUsed);
310 if((fVersionVdriftUsed == 0) && (fVersionGainUsed == 0)) fStatus = -1;
315 //___________________________________________________________________________________________________________________
317 Bool_t AliTRDPreprocessorOffline::ReadGainGlobal(const Char_t* fileName){
319 // read calibration entries from file
321 if(fCH2d) return kTRUE;
322 TFile fcalib(fileName);
323 TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
325 TH2I *ch2d = (TH2I *) array->FindObject("CH2d");
326 if(!ch2d) return kFALSE;
327 fCH2d = (TH2I*)ch2d->Clone();
328 //fNEvents = (TH1I *) array->FindObject("NEvents");
329 //fAbsoluteGain = (TH2F *) array->FindObject("AbsoluteGain");
331 TH2I *ch2d = (TH2I *) fcalib.Get("CH2d");
332 if(!ch2d) return kFALSE;
333 fCH2d = (TH2I*)ch2d->Clone();
334 //fNEvents = (TH1I *) fcalib.Get("NEvents");
335 //fAbsoluteGain = (TH2F *) fcalib.Get("AbsoluteGain");
337 fCH2d->SetDirectory(0);
338 //printf("title of CH2d %s\n",fCH2d->GetTitle());
343 //_________________________________________________________________________________________________________________
345 Bool_t AliTRDPreprocessorOffline::ReadVdriftT0Global(const Char_t* fileName){
347 // read calibration entries from file
349 if(fPH2d) return kTRUE;
350 TFile fcalib(fileName);
351 TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
353 TProfile2D *ph2d = (TProfile2D *) array->FindObject("PH2d");
354 if(!ph2d) return kFALSE;
355 fPH2d = (TProfile2D*)ph2d->Clone();
356 //fNEvents = (TH1I *) array->FindObject("NEvents");
358 TProfile2D *ph2d = (TProfile2D *) fcalib.Get("PH2d");
359 if(!ph2d) return kFALSE;
360 fPH2d = (TProfile2D*)ph2d->Clone();
361 //fNEvents = (TH1I *) fcalib.Get("NEvents");
363 fPH2d->SetDirectory(0);
364 //printf("title of PH2d %s\n",fPH2d->GetTitle());
369 //___________________________________________________________________________________________________________________
371 Bool_t AliTRDPreprocessorOffline::ReadVdriftLinearFitGlobal(const Char_t* fileName){
373 // read calibration entries from file
375 if(fAliTRDCalibraVdriftLinearFit) return kTRUE;
376 TFile fcalib(fileName);
377 TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
379 fAliTRDCalibraVdriftLinearFit = (AliTRDCalibraVdriftLinearFit *) array->FindObject("AliTRDCalibraVdriftLinearFit");
380 //fNEvents = (TH1I *) array->FindObject("NEvents");
382 fAliTRDCalibraVdriftLinearFit = (AliTRDCalibraVdriftLinearFit *) fcalib.Get("AliTRDCalibraVdriftLinearFit");
383 //fNEvents = (TH1I *) fcalib.Get("NEvents");
385 if(!fAliTRDCalibraVdriftLinearFit) {
386 //printf("No AliTRDCalibraVdriftLinearFit\n");
392 //_____________________________________________________________________________________________________________
394 Bool_t AliTRDPreprocessorOffline::ReadPRFGlobal(const Char_t* fileName){
396 // read calibration entries from file
398 if(fPRF2d) return kTRUE;
399 TFile fcalib(fileName);
400 TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
402 TProfile2D *prf2d = (TProfile2D *) array->FindObject("PRF2d");
403 if(!prf2d) return kFALSE;
404 fPRF2d = (TProfile2D*)prf2d->Clone();
405 //fNEvents = (TH1I *) array->FindObject("NEvents");
407 TProfile2D *prf2d = (TProfile2D *) fcalib.Get("PRF2d");
408 if(!prf2d) return kFALSE;
409 fPRF2d = (TProfile2D*)prf2d->Clone();
410 //fNEvents = (TH1I *) fcalib.Get("NEvents");
412 fPRF2d->SetDirectory(0);
413 //printf("title of PRF2d %s\n",fPRF2d->GetTitle());
418 //__________________________________________________________________________________________________________
420 Bool_t AliTRDPreprocessorOffline::AnalyzeGain(){
422 // Analyze gain - produce the calibration objects
425 AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
426 calibra->SetMinEntries(fMinStatsGain); // If there is less than 1000 entries in the histo: no fit
427 calibra->AnalyseCH(fCH2d);
429 Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(0))
430 + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(0));
431 Int_t nbfit = calibra->GetNumberFit();
432 Int_t nbE = calibra->GetNumberEnt();
436 Bool_t meanother = kFALSE;
439 (nbfit >= 0.5*nbE) && (nbE > 30)) {
440 // create the cal objects
441 if(!fBackCorrectGain) {
442 calibra->PutMeanValueOtherVectorFit(1,kTRUE);
445 TObjArray object = calibra->GetVectorFit();
446 AliTRDCalDet *calDetGain = calibra->CreateDetObjectGain(&object,meanother);
447 TH1F *coefGain = calDetGain->MakeHisto1DAsFunctionOfDet();
448 // Put them in the array
449 fCalibObjects->AddAt(calDetGain,kGain);
450 fPlots->AddAt(coefGain,kGain);
455 fNotEnoughStatisticsForTheGain = kTRUE;
456 Int_t minStatsGain = fMinStatsGain*30;
457 calibra->SetMinEntries(minStatsGain); // Because we do it for all, we increase this
458 Double_t gainoverallnotnormalized = calibra->AnalyseCHAllTogether(fCH2d);
459 if(fCalDetGainUsed && (gainoverallnotnormalized > 0.0)) {
460 AliTRDCalDet *calDetGain = new AliTRDCalDet(*fCalDetGainUsed);
461 Double_t oldmean = fCalDetGainUsed->CalcMean(kFALSE);
462 //printf("oldmean %f\n",oldmean);
464 Double_t scalefactor = calibra->GetScaleFactorGain();
465 //printf("Correction factor %f\n",gainoverallnotnormalized*scalefactor);
466 calDetGain->Multiply(gainoverallnotnormalized*scalefactor/oldmean);
467 //printf("newmean %f\n",calDetGain->CalcMean(kFALSE));
468 TH1F *coefGain = calDetGain->MakeHisto1DAsFunctionOfDet();
469 fCalibObjects->AddAt(calDetGain,kGain);
470 fPlots->AddAt(coefGain,kGain);
484 calibra->ResetVectorFit();
489 //_____________________________________________________________________________________________________
490 Bool_t AliTRDPreprocessorOffline::AnalyzeVdriftT0(){
492 // Analyze VdriftT0 - produce the calibration objects
495 AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
496 calibra->SetMinEntries(fMinStatsVdriftT0PH); // If there is less than 1000 entries in the histo: no fit
497 calibra->AnalysePH(fPH2d);
499 Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(1))
500 + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(1));
501 Int_t nbfit = calibra->GetNumberFit();
502 Int_t nbfitSuccess = calibra->GetNumberFitSuccess();
503 Int_t nbE = calibra->GetNumberEnt();
505 //printf("nbtg %d, nbfit %d, nbE %d, nbfitSuccess %d\n",nbtg,nbfit,nbE,nbfitSuccess);
509 (nbfit >= 0.5*nbE) && (nbE > 30) && (nbfitSuccess > 30)) {
510 //printf("Pass the cut for VdriftT0\n");
511 // create the cal objects
512 calibra->RemoveOutliers(1,kFALSE);
513 calibra->PutMeanValueOtherVectorFit(1,kFALSE);
514 calibra->RemoveOutliers2(kFALSE);
515 calibra->PutMeanValueOtherVectorFit2(1,kFALSE);
517 TObjArray object = calibra->GetVectorFit();
518 AliTRDCalDet *calDetVdrift = calibra->CreateDetObjectVdrift(&object);
519 TH1F *coefVdriftPH = calDetVdrift->MakeHisto1DAsFunctionOfDet();
520 AliTRDCalPad *calPadVdrift = (AliTRDCalPad *)calibra->CreatePadObjectVdrift(&object,calDetVdrift);
521 TH1F *coefPadVdrift = calPadVdrift->MakeHisto1D();
522 object = calibra->GetVectorFit2();
523 AliTRDCalDet *calDetT0 = calibra->CreateDetObjectT0(&object);
524 TH1F *coefT0 = calDetT0->MakeHisto1DAsFunctionOfDet();
525 AliTRDCalPad *calPadT0 = (AliTRDCalPad *)calibra->CreatePadObjectT0(&object,calDetT0);
526 TH1F *coefPadT0 = calPadT0->MakeHisto1D();
527 // Put them in the array
528 fCalibObjects->AddAt(calDetT0,kT0PHDet);
529 fCalibObjects->AddAt(calDetVdrift,kVdriftPHDet);
530 fCalibObjects->AddAt(calPadT0,kT0PHPad);
531 fCalibObjects->AddAt(calPadVdrift,kVdriftPHPad);
532 fPlots->AddAt(coefVdriftPH,kVdriftPHDet);
533 fPlots->AddAt(coefT0,kT0PHDet);
534 fPlots->AddAt(coefPadVdrift,kVdriftPHPad);
535 fPlots->AddAt(coefPadT0,kT0PHPad);
542 calibra->ResetVectorFit();
547 //____________________________________________________________________________________________________________________
548 Bool_t AliTRDPreprocessorOffline::AnalyzeVdriftLinearFit(){
550 // Analyze vdrift linear fit - produce the calibration objects
553 //printf("Analyse linear fit\n");
556 AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
557 calibra->SetMinEntries(fMinStatsVdriftLinear); // If there is less than 1000 entries in the histo: no fit
558 //printf("Fill PE Array\n");
559 fAliTRDCalibraVdriftLinearFit->FillPEArray();
560 //printf("AliTRDCalibraFit\n");
561 calibra->AnalyseLinearFitters(fAliTRDCalibraVdriftLinearFit);
565 Int_t nbfit = calibra->GetNumberFit();
566 Int_t nbE = calibra->GetNumberEnt();
571 if ((nbfit >= 0.5*nbE) && (nbE > 30)) {
572 // create the cal objects
573 //calibra->RemoveOutliers(1,kTRUE);
574 calibra->PutMeanValueOtherVectorFit(1,kTRUE);
575 //calibra->RemoveOutliers2(kTRUE);
576 calibra->PutMeanValueOtherVectorFit2(1,kTRUE);
578 TObjArray object = calibra->GetVectorFit();
579 AliTRDCalDet *calDetVdrift = calibra->CreateDetObjectVdrift(&object,kTRUE);
580 TH1F *coefDriftLinear = calDetVdrift->MakeHisto1DAsFunctionOfDet();
581 object = calibra->GetVectorFit2();
582 AliTRDCalDet *calDetLorentz = calibra->CreateDetObjectLorentzAngle(&object);
583 TH1F *coefLorentzAngle = calDetLorentz->MakeHisto1DAsFunctionOfDet();
584 //if(!calDetLorentz) printf("No lorentz created\n");
585 // Put them in the array
586 fCalibObjects->AddAt(calDetVdrift,kVdriftLinear);
587 fCalibObjects->AddAt(calDetLorentz,kLorentzLinear);
588 fPlots->AddAt(coefDriftLinear,kVdriftLinear);
589 fPlots->AddAt(coefLorentzAngle,kLorentzLinear);
594 fNotEnoughStatisticsForTheVdriftLinear = kTRUE;
595 Int_t minNumberOfEntriesForAll = fMinStatsVdriftLinear*30;
596 calibra->SetMinEntries(minNumberOfEntriesForAll); // Because we do it for all, we increase this
597 Double_t vdriftoverall = calibra->AnalyseLinearFittersAllTogether(fAliTRDCalibraVdriftLinearFit);
598 if(fCalDetVdriftUsed && (vdriftoverall > 0.0)) {
599 AliTRDCalDet *calDetVdrift = new AliTRDCalDet(*fCalDetVdriftUsed);
600 Double_t oldmean = fCalDetVdriftUsed->CalcMean(kFALSE);
601 //printf("oldmean %f\n",oldmean);
603 //printf("Correction factor %f\n",vdriftoverall);
604 calDetVdrift->Multiply(vdriftoverall/oldmean);
605 //printf("newmean %f\n",calDetVdrift->CalcMean(kFALSE));
606 TH1F *coefDriftLinear = calDetVdrift->MakeHisto1DAsFunctionOfDet();
607 // Put them in the array
608 fCalibObjects->AddAt(calDetVdrift,kVdriftLinear);
609 fPlots->AddAt(coefDriftLinear,kVdriftLinear);
619 calibra->ResetVectorFit();
624 //________________________________________________________________________________________________________________
626 Bool_t AliTRDPreprocessorOffline::AnalyzePRF(){
628 // Analyze PRF - produce the calibration objects
631 AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
632 calibra->SetMinEntries(fMinStatsPRF); // If there is less than 1000 entries in the histo: no fit
633 calibra->AnalysePRFMarianFit(fPRF2d);
635 Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(2))
636 + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(2));
637 Int_t nbfit = calibra->GetNumberFit();
638 Int_t nbE = calibra->GetNumberEnt();
644 (nbfit >= 0.95*nbE) && (nbE > 30)) {
645 // create the cal objects
646 TObjArray object = calibra->GetVectorFit();
647 AliTRDCalPad *calPadPRF = (AliTRDCalPad*) calibra->CreatePadObjectPRF(&object);
648 TH1F *coefPRF = calPadPRF->MakeHisto1D();
649 // Put them in the array
650 fCalibObjects->AddAt(calPadPRF,kPRF);
651 fPlots->AddAt(coefPRF,kPRF);
656 calibra->ResetVectorFit();
662 //_____________________________________________________________________________
663 Bool_t AliTRDPreprocessorOffline::AnalyzeChamberStatus()
666 // Produce AliTRDCalChamberStatus out of calibration results
669 // set up AliTRDCalChamberStatus
670 AliTRDCalChamberStatus *CalChamberStatus = new AliTRDCalChamberStatus();
671 for(Int_t det = 0; det < 540; det++) CalChamberStatus->SetStatus(det,1);
673 // get calibration objects
674 AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain);
675 AliTRDCalDet *calDetVDrift = (AliTRDCalDet *) fCalibObjects->At(kVdriftLinear);
678 if((!calDetGain) || (!calDetVDrift) || (!fCH2d)) return kFALSE;
681 Double_t gainmean = calDetGain->GetMean();
682 Double_t vdriftmean = calDetVDrift->GetMean();
684 Double_t gainrms = calDetGain->GetRMSRobust();
685 Double_t vdriftrms = calDetVDrift->GetRMSRobust();
687 //printf("Gain mean: %f, rms: %f\n",gainmean,gainrms);
688 //printf("Vdrift mean: %f, rms: %f\n",vdriftmean,vdriftrms);
691 if((TMath::Abs(gainrms) < 0.001) || (TMath::Abs(vdriftrms) < 0.001)) return kFALSE;
693 // mask chambers with empty gain entries
695 for (Int_t idet = 0; idet < 540; idet++) {
698 TH1I *projch = (TH1I *) fCH2d->ProjectionX("projch",idet+1,idet+1,(Option_t *)"e");
699 Double_t entries = projch->GetEntries();
702 Double_t gain = calDetGain->GetValue(idet);
705 Double_t vdrift = calDetVDrift->GetValue(idet);
709 TMath::Abs(gainmean-gain) > (15.0*gainrms) ||
710 TMath::Abs(vdriftmean-vdrift) > (15.0*vdriftrms)) {
712 //printf(" chamber det %03d masked \n",idet);
713 //printf(" gainmean %f and gain %f, gainrms %f \n",gainmean,gain,gainrms);
714 //printf(" vdriftmean %f and vdrift %f, vdriftrms %f \n",vdriftmean,vdrift,vdriftrms);
715 CalChamberStatus->SetStatus(idet,AliTRDCalChamberStatus::kMasked);
720 // installed supermodules+1 -> abort
721 if(counter > (7+1)*30) {
722 printf("ERROR: more than one SM to be masked!! \n Abort...\n");
723 if(projch) delete projch;
733 for(Int_t sm=0; sm < 18; sm++) {
735 for(Int_t det = 0; det < 30; det++){
736 Int_t detector = sm*30+det;
737 if(CalChamberStatus->IsMasked(detector)) counter++;
740 for(Int_t det = 0; det < 30; det++){
741 Int_t detector = sm*30+det;
742 CalChamberStatus->SetStatus(detector,AliTRDCalChamberStatus::kInstalled);
747 fCalibObjects->AddAt(CalChamberStatus,kChamberStatus);
753 //________________________________________________________________________________________________
754 void AliTRDPreprocessorOffline::CorrectFromDetGainUsed() {
756 // Correct from the gas gain used afterwards
758 AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain);
759 if(!calDetGain) return;
765 for(Int_t det = 0; det < 540; det++) {
767 Float_t gaininit = fCalDetGainUsed->GetValue(det);
768 Float_t gainout = calDetGain->GetValue(det);
771 if(TMath::Abs(gainout-1.0) > 0.000001) {
772 mean += (gaininit*gainout);
776 if(nbdet > 0) mean = mean/nbdet;
778 for(Int_t det = 0; det < 540; det++) {
780 Float_t gaininit = fCalDetGainUsed->GetValue(det);
781 Float_t gainout = calDetGain->GetValue(det);
783 if(TMath::Abs(gainout-1.0) > 0.000001) {
784 Double_t newgain = gaininit*gainout;
785 if(newgain < 0.1) newgain = 0.1;
786 if(newgain > 1.9) newgain = 1.9;
787 calDetGain->SetValue(det,newgain);
790 Double_t newgain = mean;
791 if(newgain < 0.1) newgain = 0.1;
792 if(newgain > 1.9) newgain = 1.9;
793 calDetGain->SetValue(det,newgain);
799 //________________________________________________________________________________________________
800 void AliTRDPreprocessorOffline::CorrectFromDetVdriftUsed() {
802 // Correct from the drift velocity
805 //printf("Correct for vdrift\n");
807 AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain);
808 if(!calDetGain) return;
810 Int_t detVdrift = kVdriftPHDet;
811 if(fMethodSecond) detVdrift = kVdriftLinear;
812 AliTRDCalDet *calDetVdrift = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
813 if(!calDetVdrift) return;
816 if(!fNotEnoughStatisticsForTheVdriftLinear) {
817 for(Int_t det = 0; det < 540; det++) {
819 Float_t vdriftinit = fCalDetVdriftUsed->GetValue(det);
820 Float_t vdriftout = calDetVdrift->GetValue(det);
822 Float_t gain = calDetGain->GetValue(det);
823 if(vdriftout > 0.0) gain = gain*vdriftinit/vdriftout;
824 if(gain < 0.1) gain = 0.1;
825 if(gain > 1.9) gain = 1.9;
826 calDetGain->SetValue(det,gain);
831 Float_t vdriftinit = fCalDetVdriftUsed->CalcMean(kFALSE);
832 Float_t vdriftout = calDetVdrift->CalcMean(kFALSE);
833 Float_t factorcorrectif = 1.0;
834 if(vdriftout > 0.0) factorcorrectif = vdriftinit/vdriftout;
835 for(Int_t det = 0; det < 540; det++) {
836 Float_t gain = calDetGain->GetValue(det);
837 gain = gain*factorcorrectif;
838 if(gain < 0.1) gain = 0.1;
839 if(gain > 1.9) gain = 1.9;
840 calDetGain->SetValue(det,gain);
846 //_________________________________________________________________________________________________________________
847 void AliTRDPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
852 AliCDBMetaData *metaData= new AliCDBMetaData();
853 metaData->SetObjectClassName("AliTRDCalDet");
854 metaData->SetResponsible("Raphaelle Bailhache");
855 metaData->SetBeamPeriod(1);
857 AliCDBId id1("TRD/Calib/ChamberGainFactor", startRunNumber, endRunNumber);
858 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
859 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain);
860 if(calDet) gStorage->Put(calDet, id1, metaData);
864 //___________________________________________________________________________________________________________________
865 void AliTRDPreprocessorOffline::UpdateOCDBExB(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
870 Int_t detExB = kLorentzLinear;
871 if(!fMethodSecond) return;
875 AliCDBMetaData *metaData= new AliCDBMetaData();
876 metaData->SetObjectClassName("AliTRDCalDet");
877 metaData->SetResponsible("Raphaelle Bailhache");
878 metaData->SetBeamPeriod(1);
880 AliCDBId id1("TRD/Calib/ChamberExB", startRunNumber, endRunNumber);
881 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
882 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detExB);
883 if(calDet) gStorage->Put(calDet, id1, metaData);
884 //if(!calDet) printf("No caldet\n");
887 //___________________________________________________________________________________________________________________
888 void AliTRDPreprocessorOffline::UpdateOCDBVdrift(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
893 Int_t detVdrift = kVdriftPHDet;
895 if(fMethodSecond) detVdrift = kVdriftLinear;
897 AliCDBMetaData *metaData= new AliCDBMetaData();
898 metaData->SetObjectClassName("AliTRDCalDet");
899 metaData->SetResponsible("Raphaelle Bailhache");
900 metaData->SetBeamPeriod(1);
902 AliCDBId id1("TRD/Calib/ChamberVdrift", startRunNumber, endRunNumber);
903 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
904 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
905 if(calDet) gStorage->Put(calDet, id1, metaData);
911 AliCDBMetaData *metaDataPad= new AliCDBMetaData();
912 metaDataPad->SetObjectClassName("AliTRDCalPad");
913 metaDataPad->SetResponsible("Raphaelle Bailhache");
914 metaDataPad->SetBeamPeriod(1);
916 AliCDBId id1Pad("TRD/Calib/LocalVdrift", startRunNumber, endRunNumber);
917 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad);
918 if(calPad) gStorage->Put(calPad, id1Pad, metaDataPad);
923 //________________________________________________________________________________________________________________________
924 void AliTRDPreprocessorOffline::UpdateOCDBT0(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
929 AliCDBMetaData *metaData= new AliCDBMetaData();
930 metaData->SetObjectClassName("AliTRDCalDet");
931 metaData->SetResponsible("Raphaelle Bailhache");
932 metaData->SetBeamPeriod(1);
934 AliCDBId id1("TRD/Calib/ChamberT0", startRunNumber, endRunNumber);
935 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
936 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet);
937 if(calDet) gStorage->Put(calDet, id1, metaData);
941 AliCDBMetaData *metaDataPad= new AliCDBMetaData();
942 metaDataPad->SetObjectClassName("AliTRDCalPad");
943 metaDataPad->SetResponsible("Raphaelle Bailhache");
944 metaDataPad->SetBeamPeriod(1);
946 AliCDBId id1Pad("TRD/Calib/LocalT0", startRunNumber, endRunNumber);
947 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad);
948 if(calPad) gStorage->Put(calPad, id1Pad, metaDataPad);
953 //_________________________________________________________________________________________________________________
954 void AliTRDPreprocessorOffline::UpdateOCDBPRF(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
959 AliCDBMetaData *metaData= new AliCDBMetaData();
960 metaData->SetObjectClassName("AliTRDCalPad");
961 metaData->SetResponsible("Raphaelle Bailhache");
962 metaData->SetBeamPeriod(1);
964 AliCDBId id1("TRD/Calib/PRFWidth", startRunNumber, endRunNumber);
965 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
966 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF);
967 if(calPad) gStorage->Put(calPad, id1, metaData);
971 //_________________________________________________________________________________________________________________
972 void AliTRDPreprocessorOffline::UpdateOCDBChamberStatus(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
977 AliCDBMetaData *metaData= new AliCDBMetaData();
978 metaData->SetObjectClassName("AliTRDCalChamberStatus");
979 metaData->SetResponsible("Raphaelle Bailhache");
980 metaData->SetBeamPeriod(1);
982 AliCDBId id1("TRD/Calib/ChamberStatus", startRunNumber, endRunNumber);
983 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
984 AliTRDCalChamberStatus *calChamberStatus = (AliTRDCalChamberStatus *) fCalibObjects->At(kChamberStatus);
985 if(calChamberStatus) gStorage->Put(calChamberStatus, id1, metaData);
989 //__________________________________________________________________________________________________________________________
990 Bool_t AliTRDPreprocessorOffline::ValidateGain() {
992 // Validate OCDB entry
995 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain);
997 Double_t mean = calDet->GetMean();
998 Double_t rms = calDet->GetRMSRobust();
999 if((mean > 0.2) && (mean < 1.4) && (rms < 0.5)) return kTRUE;
1000 //if((mean > 0.2) && (mean < 1.4)) return kTRUE;
1002 fStatus += 1000000000;
1011 //__________________________________________________________________________________________________________________________
1012 Bool_t AliTRDPreprocessorOffline::ValidateVdrift(){
1014 // Update OCDB entry
1017 Int_t detVdrift = kVdriftPHDet;
1020 if(fMethodSecond) detVdrift = kVdriftLinear;
1022 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
1024 Double_t mean = calDet->GetMean();
1025 Double_t rms = calDet->GetRMSRobust();
1026 //printf("Vdrift::mean %f, rms %f\n",mean,rms);
1027 if(!((mean > 1.0) && (mean < 2.0) && (rms < 0.5))) {
1028 fStatus += 10000000;
1034 if(!fMethodSecond) {
1035 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad);
1037 Double_t mean = calPad->GetMean();
1038 Double_t rms = calPad->GetRMS();
1039 //printf("Vdrift::meanpad %f, rmspad %f\n",mean,rms);
1040 if(!((mean > 0.9) && (mean < 1.1) && (rms < 0.6))) {
1041 fStatus += 10000000;
1051 //__________________________________________________________________________________________________________________________
1052 Bool_t AliTRDPreprocessorOffline::ValidateT0(){
1054 // Update OCDB entry
1057 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet);
1058 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad);
1059 if(calDet && calPad) {
1060 Double_t meandet = calDet->GetMean();
1061 Double_t rmsdet = calDet->GetRMSRobust();
1062 Double_t meanpad = calPad->GetMean();
1063 //Double_t rmspad = calPad->GetRMS();
1064 //printf("T0::minimum %f, rmsdet %f,meanpad %f, rmspad %f\n",meandet,rmsdet,meanpad,rmspad);
1065 if((meandet > -1.5) && (meandet < 5.0) && (rmsdet < 4.0) && (meanpad < 5.0) && (meanpad > -0.5)) return kTRUE;
1067 fStatus += 100000000;
1074 //__________________________________________________________________________________________________________________________
1075 Bool_t AliTRDPreprocessorOffline::ValidatePRF() const{
1077 // Update OCDB entry
1080 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF);
1082 Double_t meanpad = calPad->GetMean();
1083 Double_t rmspad = calPad->GetRMS();
1084 //printf("PRF::meanpad %f, rmspad %f\n",meanpad,rmspad);
1085 if((meanpad < 1.0) && (rmspad < 0.8)) return kTRUE;
1092 //__________________________________________________________________________________________________________________________
1093 Bool_t AliTRDPreprocessorOffline::ValidateChamberStatus() const{
1095 // Update OCDB entry
1098 AliTRDCalChamberStatus *calChamberStatus = (AliTRDCalChamberStatus *) fCalibObjects->At(kChamberStatus);
1099 if(calChamberStatus) {
1100 Int_t detectormasked = 0;
1101 for(Int_t det = 0; det < 540; det++) {
1102 if(calChamberStatus->IsMasked(det)) detectormasked++;
1104 //printf("Number of chambers masked %d\n",detectormasked);
1105 if(detectormasked > 40) return kFALSE;
1111 //_____________________________________________________________________________
1112 Int_t AliTRDPreprocessorOffline::GetVersion(TString name) const
1115 // Get version from the title
1119 const Char_t *version = "Ver";
1120 if(!strstr(name.Data(),version)) return -1;
1121 const Char_t *after = "Subver";
1122 if(!strstr(name.Data(),after)) return -1;
1124 for(Int_t ver = 0; ver < 999999999; ver++) {
1126 TString vertry(version);
1130 //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
1132 if(strstr(name.Data(),vertry.Data())) return ver;
1140 //_____________________________________________________________________________
1141 Int_t AliTRDPreprocessorOffline::GetSubVersion(TString name) const
1144 // Get subversion from the title
1148 const Char_t *subversion = "Subver";
1149 if(!strstr(name.Data(),subversion)) return -1;
1150 const Char_t *after = "FirstRun";
1151 if(!strstr(name.Data(),after)) {
1154 if(!strstr(name.Data(),after)) return -1;
1157 for(Int_t ver = 0; ver < 999999999; ver++) {
1159 TString vertry(subversion);
1163 //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
1165 if(strstr(name.Data(),vertry.Data())) return ver;
1173 //_____________________________________________________________________________
1174 Int_t AliTRDPreprocessorOffline::GetFirstRun(TString name) const
1177 // Get first run from the title
1181 const Char_t *firstrun = "FirstRun";
1182 if(!strstr(name.Data(),firstrun)) return -1;
1183 const Char_t *after = "Nz";
1184 if(!strstr(name.Data(),after)) return -1;
1187 for(Int_t ver = 0; ver < 999999999; ver++) {
1189 TString vertry(firstrun);
1193 //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
1195 if(strstr(name.Data(),vertry.Data())) return ver;