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"
62 ClassImp(AliTRDPreprocessorOffline)
64 AliTRDPreprocessorOffline::AliTRDPreprocessorOffline():
65 TNamed("TPCPreprocessorOffline","TPCPreprocessorOffline"),
67 fNameList("TRDCalib"),
69 fCalDetVdriftUsed(0x0),
73 fAliTRDCalibraVdriftLinearFit(0x0),
76 fPlots(new TObjArray(8)),
77 fCalibObjects(new TObjArray(8)),
79 fSubVersionGainUsed(0),
80 fVersionVdriftUsed(0),
81 fSubVersionVdriftUsed(0),
82 fSwitchOnValidation(kTRUE),
83 fVdriftValidated(kFALSE),
85 fMinStatsVdriftT0PH(800*20),
86 fMinStatsVdriftLinear(800),
91 // default constructor
94 //_________________________________________________________________________________________________________________
95 AliTRDPreprocessorOffline::~AliTRDPreprocessorOffline() {
100 if(fCalDetGainUsed) delete fCalDetGainUsed;
101 if(fCalDetVdriftUsed) delete fCalDetVdriftUsed;
102 if(fCH2d) delete fCH2d;
103 if(fPH2d) delete fPH2d;
104 if(fPRF2d) delete fPRF2d;
105 if(fAliTRDCalibraVdriftLinearFit) delete fAliTRDCalibraVdriftLinearFit;
106 if(fNEvents) delete fNEvents;
107 if(fAbsoluteGain) delete fAbsoluteGain;
108 if(fPlots) delete fPlots;
109 if(fCalibObjects) delete fCalibObjects;
112 //___________________________________________________________________________________________________________________
114 void AliTRDPreprocessorOffline::CalibVdriftT0(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){
116 // make calibration of the drift velocity
118 // file - the location of input file
119 // startRunNumber, endRunNumber - run validity period
120 // ocdbStorage - path to the OCDB storage
121 // - if empty - local storage 'pwd' uesed
122 if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
126 fVdriftValidated = kTRUE;
127 fT0Validated = kTRUE;
129 // 2. extraction of the information
131 if(ReadVdriftT0Global(file)) AnalyzeVdriftT0();
132 if(ReadVdriftLinearFitGlobal(file)) AnalyzeVdriftLinearFit();
134 // 3. Append QA plots
136 //MakeDefaultPlots(fVdriftArray,fVdriftArray);
139 // 4. validate OCDB entries
141 if(fSwitchOnValidation==kTRUE && ValidateVdrift()==kFALSE) {
142 AliError("TRD vdrift OCDB parameters out of range!");
143 fVdriftValidated = kFALSE;
145 if(fSwitchOnValidation==kTRUE && ValidateT0()==kFALSE) {
146 AliError("TRD t0 OCDB parameters out of range!");
147 fT0Validated = kFALSE;
153 if(fVdriftValidated) UpdateOCDBVdrift(startRunNumber,endRunNumber,ocdbStorage);
154 if(fT0Validated) UpdateOCDBT0(startRunNumber,endRunNumber,ocdbStorage);
157 //_________________________________________________________________________________________________________________
159 void AliTRDPreprocessorOffline::CalibGain(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){
161 // make calibration of the drift velocity
163 // file - the location of input file
164 // startRunNumber, endRunNumber - run validity period
165 // ocdbStorage - path to the OCDB storage
166 // - if empty - local storage 'pwd' uesed
167 if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
170 if(!ReadGainGlobal(file)) return;
173 // 2. extraction of the information
176 if(fCalDetGainUsed) CorrectFromDetGainUsed();
177 if(fCalDetVdriftUsed) CorrectFromDetVdriftUsed();
179 // 3. Append QA plots
181 //MakeDefaultPlots(fVdriftArray,fVdriftArray);
184 // 4. validate OCDB entries
186 if(fSwitchOnValidation==kTRUE && ValidateGain()==kFALSE) {
187 AliError("TRD gain OCDB parameters out of range!");
194 if((!fCalDetVdriftUsed) || (fCalDetVdriftUsed && fVdriftValidated)) UpdateOCDBGain(startRunNumber,endRunNumber,ocdbStorage);
198 //________________________________________________________________________________________________________________
200 void AliTRDPreprocessorOffline::CalibPRF(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){
202 // make calibration of the drift velocity
204 // file - the location of input file
205 // startRunNumber, endRunNumber - run validity period
206 // ocdbStorage - path to the OCDB storage
207 // - if empty - local storage 'pwd' uesed
208 if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
211 if(!ReadPRFGlobal(file)) return;
214 // 2. extraction of the information
218 // 3. Append QA plots
220 //MakeDefaultPlots(fVdriftArray,fVdriftArray);
224 // 4. validate OCDB entries
226 if(fSwitchOnValidation==kTRUE && ValidatePRF()==kFALSE) {
227 AliError("TRD prf OCDB parameters out of range!");
234 UpdateOCDBPRF(startRunNumber,endRunNumber,ocdbStorage);
237 //______________________________________________________________________________________________________
238 Bool_t AliTRDPreprocessorOffline::Init(const Char_t* fileName){
240 // read the calibration used during the reconstruction
243 if(ReadVdriftT0Global(fileName)) {
245 TString nameph = fPH2d->GetTitle();
246 fVersionVdriftUsed = GetVersion(nameph);
247 fSubVersionVdriftUsed = GetSubVersion(nameph);
249 //printf("Found Version %d, Subversion %d for vdrift\n",fVersionVdriftUsed,fSubVersionVdriftUsed);
253 if(ReadGainGlobal(fileName)) {
255 TString namech = fCH2d->GetTitle();
256 fVersionGainUsed = GetVersion(namech);
257 fSubVersionGainUsed = GetSubVersion(namech);
259 //printf("Found Version %d, Subversion %d for gain\n",fVersionGainUsed,fSubVersionGainUsed);
266 //___________________________________________________________________________________________________________________
268 Bool_t AliTRDPreprocessorOffline::ReadGainGlobal(const Char_t* fileName){
270 // read calibration entries from file
272 if(fCH2d) return kTRUE;
273 TFile fcalib(fileName);
274 TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
276 TH2I *ch2d = (TH2I *) array->FindObject("CH2d");
277 if(!ch2d) return kFALSE;
278 fCH2d = (TH2I*)ch2d->Clone();
279 //fNEvents = (TH1I *) array->FindObject("NEvents");
280 //fAbsoluteGain = (TH2F *) array->FindObject("AbsoluteGain");
282 TH2I *ch2d = (TH2I *) fcalib.Get("CH2d");
283 if(!ch2d) return kFALSE;
284 fCH2d = (TH2I*)ch2d->Clone();
285 //fNEvents = (TH1I *) fcalib.Get("NEvents");
286 //fAbsoluteGain = (TH2F *) fcalib.Get("AbsoluteGain");
288 fCH2d->SetDirectory(0);
289 //printf("title of CH2d %s\n",fCH2d->GetTitle());
294 //_________________________________________________________________________________________________________________
296 Bool_t AliTRDPreprocessorOffline::ReadVdriftT0Global(const Char_t* fileName){
298 // read calibration entries from file
300 if(fPH2d) return kTRUE;
301 TFile fcalib(fileName);
302 TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
304 TProfile2D *ph2d = (TProfile2D *) array->FindObject("PH2d");
305 if(!ph2d) return kFALSE;
306 fPH2d = (TProfile2D*)ph2d->Clone();
307 //fNEvents = (TH1I *) array->FindObject("NEvents");
309 TProfile2D *ph2d = (TProfile2D *) fcalib.Get("PH2d");
310 if(!ph2d) return kFALSE;
311 fPH2d = (TProfile2D*)ph2d->Clone();
312 //fNEvents = (TH1I *) fcalib.Get("NEvents");
314 fPH2d->SetDirectory(0);
315 //printf("title of PH2d %s\n",fPH2d->GetTitle());
320 //___________________________________________________________________________________________________________________
322 Bool_t AliTRDPreprocessorOffline::ReadVdriftLinearFitGlobal(const Char_t* fileName){
324 // read calibration entries from file
326 if(fAliTRDCalibraVdriftLinearFit) return kTRUE;
327 TFile fcalib(fileName);
328 TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
330 fAliTRDCalibraVdriftLinearFit = (AliTRDCalibraVdriftLinearFit *) array->FindObject("AliTRDCalibraVdriftLinearFit");
331 //fNEvents = (TH1I *) array->FindObject("NEvents");
333 fAliTRDCalibraVdriftLinearFit = (AliTRDCalibraVdriftLinearFit *) fcalib.Get("AliTRDCalibraVdriftLinearFit");
334 //fNEvents = (TH1I *) fcalib.Get("NEvents");
336 if(!fAliTRDCalibraVdriftLinearFit) {
337 //printf("No AliTRDCalibraVdriftLinearFit\n");
343 //_____________________________________________________________________________________________________________
345 Bool_t AliTRDPreprocessorOffline::ReadPRFGlobal(const Char_t* fileName){
347 // read calibration entries from file
349 if(fPRF2d) return kTRUE;
350 TFile fcalib(fileName);
351 TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
353 TProfile2D *prf2d = (TProfile2D *) array->FindObject("PRF2d");
354 if(!prf2d) return kFALSE;
355 fPRF2d = (TProfile2D*)prf2d->Clone();
356 //fNEvents = (TH1I *) array->FindObject("NEvents");
358 TProfile2D *prf2d = (TProfile2D *) fcalib.Get("PRF2d");
359 if(!prf2d) return kFALSE;
360 fPRF2d = (TProfile2D*)prf2d->Clone();
361 //fNEvents = (TH1I *) fcalib.Get("NEvents");
363 fPRF2d->SetDirectory(0);
364 //printf("title of PRF2d %s\n",fPRF2d->GetTitle());
369 //__________________________________________________________________________________________________________
371 Bool_t AliTRDPreprocessorOffline::AnalyzeGain(){
373 // Analyze gain - produce the calibration objects
376 AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
377 calibra->SetMinEntries(fMinStatsGain); // If there is less than 1000 entries in the histo: no fit
378 calibra->AnalyseCH(fCH2d);
380 Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(0))
381 + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(0));
382 Int_t nbfit = calibra->GetNumberFit();
383 Int_t nbE = calibra->GetNumberEnt();
387 Bool_t meanother = kFALSE;
390 (nbfit >= 0.5*nbE) && (nbE > 30)) {
391 // create the cal objects
392 if(!fCalDetGainUsed) {
393 calibra->PutMeanValueOtherVectorFit(1,kTRUE);
396 TObjArray object = calibra->GetVectorFit();
397 AliTRDCalDet *calDetGain = calibra->CreateDetObjectGain(&object,meanother);
398 TH1F *coefGain = calDetGain->MakeHisto1DAsFunctionOfDet();
399 // Put them in the array
400 fCalibObjects->AddAt(calDetGain,kGain);
401 fPlots->AddAt(coefGain,kGain);
406 calibra->ResetVectorFit();
411 //_____________________________________________________________________________________________________
412 Bool_t AliTRDPreprocessorOffline::AnalyzeVdriftT0(){
414 // Analyze VdriftT0 - produce the calibration objects
417 AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
418 calibra->SetMinEntries(fMinStatsVdriftT0PH); // If there is less than 1000 entries in the histo: no fit
419 calibra->AnalysePH(fPH2d);
421 Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(1))
422 + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(1));
423 Int_t nbfit = calibra->GetNumberFit();
424 Int_t nbfitSuccess = calibra->GetNumberFitSuccess();
425 Int_t nbE = calibra->GetNumberEnt();
427 //printf("nbtg %d, nbfit %d, nbE %d, nbfitSuccess %d\n",nbtg,nbfit,nbE,nbfitSuccess);
431 (nbfit >= 0.5*nbE) && (nbE > 30) && (nbfitSuccess > 30)) {
432 //printf("Pass the cut for VdriftT0\n");
433 // create the cal objects
434 calibra->RemoveOutliers(1,kFALSE);
435 calibra->PutMeanValueOtherVectorFit(1,kFALSE);
436 calibra->RemoveOutliers2(kFALSE);
437 calibra->PutMeanValueOtherVectorFit2(1,kFALSE);
439 TObjArray object = calibra->GetVectorFit();
440 AliTRDCalDet *calDetVdrift = calibra->CreateDetObjectVdrift(&object);
441 TH1F *coefVdriftPH = calDetVdrift->MakeHisto1DAsFunctionOfDet();
442 AliTRDCalPad *calPadVdrift = (AliTRDCalPad *)calibra->CreatePadObjectVdrift(&object,calDetVdrift);
443 TH1F *coefPadVdrift = calPadVdrift->MakeHisto1D();
444 object = calibra->GetVectorFit2();
445 AliTRDCalDet *calDetT0 = calibra->CreateDetObjectT0(&object);
446 TH1F *coefT0 = calDetT0->MakeHisto1DAsFunctionOfDet();
447 AliTRDCalPad *calPadT0 = (AliTRDCalPad *)calibra->CreatePadObjectT0(&object,calDetT0);
448 TH1F *coefPadT0 = calPadT0->MakeHisto1D();
449 // Put them in the array
450 fCalibObjects->AddAt(calDetT0,kT0PHDet);
451 fCalibObjects->AddAt(calDetVdrift,kVdriftPHDet);
452 fCalibObjects->AddAt(calPadT0,kT0PHPad);
453 fCalibObjects->AddAt(calPadVdrift,kVdriftPHPad);
454 fPlots->AddAt(coefVdriftPH,kVdriftPHDet);
455 fPlots->AddAt(coefT0,kT0PHDet);
456 fPlots->AddAt(coefPadVdrift,kVdriftPHPad);
457 fPlots->AddAt(coefPadT0,kT0PHPad);
461 calibra->ResetVectorFit();
466 //____________________________________________________________________________________________________________________
467 Bool_t AliTRDPreprocessorOffline::AnalyzeVdriftLinearFit(){
469 // Analyze vdrift linear fit - produce the calibration objects
472 AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
473 calibra->SetMinEntries(fMinStatsVdriftLinear); // If there is less than 1000 entries in the histo: no fit
474 fAliTRDCalibraVdriftLinearFit->FillPEArray();
475 calibra->AnalyseLinearFitters(fAliTRDCalibraVdriftLinearFit);
478 Int_t nbfit = calibra->GetNumberFit();
479 Int_t nbE = calibra->GetNumberEnt();
484 if ((nbfit >= 0.5*nbE) && (nbE > 30)) {
485 // create the cal objects
486 //calibra->RemoveOutliers(1,kTRUE);
487 calibra->PutMeanValueOtherVectorFit(1,kTRUE);
488 //calibra->RemoveOutliers2(kTRUE);
489 calibra->PutMeanValueOtherVectorFit2(1,kTRUE);
491 TObjArray object = calibra->GetVectorFit();
492 AliTRDCalDet *calDetVdrift = calibra->CreateDetObjectVdrift(&object,kTRUE);
493 TH1F *coefDriftLinear = calDetVdrift->MakeHisto1DAsFunctionOfDet();
494 object = calibra->GetVectorFit2();
495 AliTRDCalDet *calDetLorentz = calibra->CreateDetObjectLorentzAngle(&object);
496 TH1F *coefLorentzAngle = calDetLorentz->MakeHisto1DAsFunctionOfDet();
497 // Put them in the array
498 fCalibObjects->AddAt(calDetVdrift,kVdriftLinear);
499 fCalibObjects->AddAt(calDetLorentz,kLorentzLinear);
500 fPlots->AddAt(coefDriftLinear,kVdriftLinear);
501 fPlots->AddAt(coefLorentzAngle,kLorentzLinear);
506 calibra->ResetVectorFit();
511 //________________________________________________________________________________________________________________
513 Bool_t AliTRDPreprocessorOffline::AnalyzePRF(){
515 // Analyze PRF - produce the calibration objects
518 AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
519 calibra->SetMinEntries(fMinStatsPRF); // If there is less than 1000 entries in the histo: no fit
520 calibra->AnalysePRFMarianFit(fPRF2d);
522 Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(2))
523 + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(2));
524 Int_t nbfit = calibra->GetNumberFit();
525 Int_t nbE = calibra->GetNumberEnt();
531 (nbfit >= 0.95*nbE) && (nbE > 30)) {
532 // create the cal objects
533 TObjArray object = calibra->GetVectorFit();
534 AliTRDCalPad *calPadPRF = (AliTRDCalPad*) calibra->CreatePadObjectPRF(&object);
535 TH1F *coefPRF = calPadPRF->MakeHisto1D();
536 // Put them in the array
537 fCalibObjects->AddAt(calPadPRF,kPRF);
538 fPlots->AddAt(coefPRF,kPRF);
543 calibra->ResetVectorFit();
548 //________________________________________________________________________________________________
549 void AliTRDPreprocessorOffline::CorrectFromDetGainUsed() {
551 // Correct from the gas gain used afterwards
553 AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain);
554 if(!calDetGain) return;
560 for(Int_t det = 0; det < 540; det++) {
562 Float_t gaininit = fCalDetGainUsed->GetValue(det);
563 Float_t gainout = calDetGain->GetValue(det);
566 if(TMath::Abs(gainout-1.0) > 0.000001) {
567 mean += (gaininit*gainout);
571 if(nbdet > 0) mean = mean/nbdet;
573 for(Int_t det = 0; det < 540; det++) {
575 Float_t gaininit = fCalDetGainUsed->GetValue(det);
576 Float_t gainout = calDetGain->GetValue(det);
578 if(TMath::Abs(gainout-1.0) > 0.000001) calDetGain->SetValue(det,gaininit*gainout);
579 else calDetGain->SetValue(det,mean);
584 //________________________________________________________________________________________________
585 void AliTRDPreprocessorOffline::CorrectFromDetVdriftUsed() {
587 // Correct from the drift velocity
590 //printf("Correct for vdrift\n");
592 AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain);
593 if(!calDetGain) return;
595 Int_t detVdrift = kVdriftPHDet;
596 if(fMethodSecond) detVdrift = kVdriftLinear;
597 AliTRDCalDet *calDetVdrift = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
598 if(!calDetVdrift) return;
602 for(Int_t det = 0; det < 540; det++) {
604 Float_t vdriftinit = fCalDetVdriftUsed->GetValue(det);
605 Float_t vdriftout = calDetVdrift->GetValue(det);
607 Float_t gain = calDetGain->GetValue(det);
608 if(vdriftout > 0.0) gain = gain*vdriftinit/vdriftout;
609 calDetGain->SetValue(det,gain);
615 //_________________________________________________________________________________________________________________
616 void AliTRDPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
621 AliCDBMetaData *metaData= new AliCDBMetaData();
622 metaData->SetObjectClassName("AliTRDCalDet");
623 metaData->SetResponsible("Raphaelle Bailhache");
624 metaData->SetBeamPeriod(1);
626 AliCDBId id1("TRD/Calib/ChamberGainFactor", startRunNumber, endRunNumber);
627 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
628 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain);
629 if(calDet) gStorage->Put(calDet, id1, metaData);
633 //___________________________________________________________________________________________________________________
634 void AliTRDPreprocessorOffline::UpdateOCDBVdrift(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
639 Int_t detVdrift = kVdriftPHDet;
641 if(fMethodSecond) detVdrift = kVdriftLinear;
643 AliCDBMetaData *metaData= new AliCDBMetaData();
644 metaData->SetObjectClassName("AliTRDCalDet");
645 metaData->SetResponsible("Raphaelle Bailhache");
646 metaData->SetBeamPeriod(1);
648 AliCDBId id1("TRD/Calib/ChamberVdrift", startRunNumber, endRunNumber);
649 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
650 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
651 if(calDet) gStorage->Put(calDet, id1, metaData);
657 AliCDBMetaData *metaDataPad= new AliCDBMetaData();
658 metaDataPad->SetObjectClassName("AliTRDCalPad");
659 metaDataPad->SetResponsible("Raphaelle Bailhache");
660 metaDataPad->SetBeamPeriod(1);
662 AliCDBId id1Pad("TRD/Calib/LocalVdrift", startRunNumber, endRunNumber);
663 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad);
664 if(calPad) gStorage->Put(calPad, id1Pad, metaDataPad);
669 //________________________________________________________________________________________________________________________
670 void AliTRDPreprocessorOffline::UpdateOCDBT0(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
675 AliCDBMetaData *metaData= new AliCDBMetaData();
676 metaData->SetObjectClassName("AliTRDCalDet");
677 metaData->SetResponsible("Raphaelle Bailhache");
678 metaData->SetBeamPeriod(1);
680 AliCDBId id1("TRD/Calib/ChamberT0", startRunNumber, endRunNumber);
681 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
682 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet);
683 if(calDet) gStorage->Put(calDet, id1, metaData);
687 AliCDBMetaData *metaDataPad= new AliCDBMetaData();
688 metaDataPad->SetObjectClassName("AliTRDCalPad");
689 metaDataPad->SetResponsible("Raphaelle Bailhache");
690 metaDataPad->SetBeamPeriod(1);
692 AliCDBId id1Pad("TRD/Calib/LocalT0", startRunNumber, endRunNumber);
693 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad);
694 if(calPad) gStorage->Put(calPad, id1Pad, metaDataPad);
699 //_________________________________________________________________________________________________________________
700 void AliTRDPreprocessorOffline::UpdateOCDBPRF(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
705 AliCDBMetaData *metaData= new AliCDBMetaData();
706 metaData->SetObjectClassName("AliTRDCalPad");
707 metaData->SetResponsible("Raphaelle Bailhache");
708 metaData->SetBeamPeriod(1);
710 AliCDBId id1("TRD/Calib/PRFWidth", startRunNumber, endRunNumber);
711 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
712 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF);
713 if(calPad) gStorage->Put(calPad, id1, metaData);
717 //__________________________________________________________________________________________________________________________
718 Bool_t AliTRDPreprocessorOffline::ValidateGain() const {
720 // Validate OCDB entry
723 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain);
725 Double_t mean = calDet->GetMean();
726 Double_t rms = calDet->GetRMS();
727 if((mean > 0.2) && (mean < 1.4) && (rms < 0.5)) return kTRUE;
734 //__________________________________________________________________________________________________________________________
735 Bool_t AliTRDPreprocessorOffline::ValidateVdrift(){
740 Int_t detVdrift = kVdriftPHDet;
743 if(fMethodSecond) detVdrift = kVdriftLinear;
745 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
747 Double_t mean = calDet->GetMean();
748 Double_t rms = calDet->GetRMS();
749 //printf("Vdrift::mean %f, rms %f\n",mean,rms);
750 if(!((mean > 1.0) && (mean < 2.0) && (rms < 0.5))) ok = kFALSE;
755 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad);
757 Double_t mean = calPad->GetMean();
758 Double_t rms = calPad->GetRMS();
759 //printf("Vdrift::meanpad %f, rmspad %f\n",mean,rms);
760 if(!((mean > 0.9) && (mean < 1.1) && (rms < 0.6))) ok = kFALSE;
768 //__________________________________________________________________________________________________________________________
769 Bool_t AliTRDPreprocessorOffline::ValidateT0(){
774 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet);
775 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad);
776 if(calDet && calPad) {
777 Double_t meandet = calDet->GetMean();
778 Double_t rmsdet = calDet->GetRMS();
779 Double_t meanpad = calPad->GetMean();
780 //Double_t rmspad = calPad->GetRMS();
781 //printf("T0::minimum %f, rmsdet %f,meanpad %f, rmspad %f\n",meandet,rmsdet,meanpad,rmspad);
782 if((meandet > -1.5) && (meandet < 5.0) && (rmsdet < 4.0) && (meanpad < 5.0) && (meanpad > -0.5)) return kTRUE;
788 //__________________________________________________________________________________________________________________________
789 Bool_t AliTRDPreprocessorOffline::ValidatePRF() const{
794 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF);
796 Double_t meanpad = calPad->GetMean();
797 Double_t rmspad = calPad->GetRMS();
798 //printf("PRF::meanpad %f, rmspad %f\n",meanpad,rmspad);
799 if((meanpad < 1.0) && (rmspad < 0.8)) return kTRUE;
806 //_____________________________________________________________________________
807 Int_t AliTRDPreprocessorOffline::GetVersion(TString name) const
810 // Get version from the title
814 const Char_t *version = "Ver";
815 if(!strstr(name.Data(),version)) return -1;
817 for(Int_t ver = 0; ver < 999999999; ver++) {
819 TString vertry(version);
823 //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
825 if(strstr(name.Data(),vertry.Data())) return ver;
833 //_____________________________________________________________________________
834 Int_t AliTRDPreprocessorOffline::GetSubVersion(TString name) const
837 // Get subversion from the title
841 const Char_t *subversion = "Subver";
842 if(!strstr(name.Data(),subversion)) return -1;
844 for(Int_t ver = 0; ver < 999999999; ver++) {
846 TString vertry(subversion);
850 //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
852 if(strstr(name.Data(),vertry.Data())) return ver;