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),
93 // default constructor
96 //_________________________________________________________________________________________________________________
97 AliTRDPreprocessorOffline::~AliTRDPreprocessorOffline() {
102 if(fCalDetGainUsed) delete fCalDetGainUsed;
103 if(fCalDetVdriftUsed) delete fCalDetVdriftUsed;
104 if(fCH2d) delete fCH2d;
105 if(fPH2d) delete fPH2d;
106 if(fPRF2d) delete fPRF2d;
107 if(fAliTRDCalibraVdriftLinearFit) delete fAliTRDCalibraVdriftLinearFit;
108 if(fNEvents) delete fNEvents;
109 if(fAbsoluteGain) delete fAbsoluteGain;
110 if(fPlots) delete fPlots;
111 if(fCalibObjects) delete fCalibObjects;
114 //___________________________________________________________________________________________________________________
116 void AliTRDPreprocessorOffline::CalibVdriftT0(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){
118 // make calibration of the drift velocity
120 // file - the location of input file
121 // startRunNumber, endRunNumber - run validity period
122 // ocdbStorage - path to the OCDB storage
123 // - if empty - local storage 'pwd' uesed
124 if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
128 fVdriftValidated = kTRUE;
129 fT0Validated = kTRUE;
131 // 2. extraction of the information
133 if(ReadVdriftT0Global(file)) AnalyzeVdriftT0();
134 //printf("Finish PH\n");
135 if(ReadVdriftLinearFitGlobal(file)) AnalyzeVdriftLinearFit();
137 // 3. Append QA plots
139 //MakeDefaultPlots(fVdriftArray,fVdriftArray);
142 // 4. validate OCDB entries
144 if(fSwitchOnValidation==kTRUE && ValidateVdrift()==kFALSE) {
145 AliError("TRD vdrift OCDB parameters out of range!");
146 fVdriftValidated = kFALSE;
148 if(fSwitchOnValidation==kTRUE && ValidateT0()==kFALSE) {
149 AliError("TRD t0 OCDB parameters out of range!");
150 fT0Validated = kFALSE;
156 if(fVdriftValidated) UpdateOCDBVdrift(startRunNumber,endRunNumber,ocdbStorage);
157 if(fT0Validated) UpdateOCDBT0(startRunNumber,endRunNumber,ocdbStorage);
160 //_________________________________________________________________________________________________________________
162 void AliTRDPreprocessorOffline::CalibGain(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){
164 // make calibration of the drift velocity
166 // file - the location of input file
167 // startRunNumber, endRunNumber - run validity period
168 // ocdbStorage - path to the OCDB storage
169 // - if empty - local storage 'pwd' uesed
170 if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
173 if(!ReadGainGlobal(file)) return;
176 // 2. extraction of the information
179 if(fCalDetGainUsed) CorrectFromDetGainUsed();
180 if(fCalDetVdriftUsed) CorrectFromDetVdriftUsed();
182 // 3. Append QA plots
184 //MakeDefaultPlots(fVdriftArray,fVdriftArray);
187 // 4. validate OCDB entries
189 if(fSwitchOnValidation==kTRUE && ValidateGain()==kFALSE) {
190 AliError("TRD gain OCDB parameters out of range!");
197 if((!fCalDetVdriftUsed) || (fCalDetVdriftUsed && fVdriftValidated)) UpdateOCDBGain(startRunNumber,endRunNumber,ocdbStorage);
201 //________________________________________________________________________________________________________________
203 void AliTRDPreprocessorOffline::CalibPRF(const Char_t* file, Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){
205 // make calibration of the drift velocity
207 // file - the location of input file
208 // startRunNumber, endRunNumber - run validity period
209 // ocdbStorage - path to the OCDB storage
210 // - if empty - local storage 'pwd' uesed
211 if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
214 if(!ReadPRFGlobal(file)) return;
217 // 2. extraction of the information
221 // 3. Append QA plots
223 //MakeDefaultPlots(fVdriftArray,fVdriftArray);
227 // 4. validate OCDB entries
229 if(fSwitchOnValidation==kTRUE && ValidatePRF()==kFALSE) {
230 AliError("TRD prf OCDB parameters out of range!");
237 UpdateOCDBPRF(startRunNumber,endRunNumber,ocdbStorage);
240 //________________________________________________________________________________________________________________
242 void AliTRDPreprocessorOffline::CalibChamberStatus(Int_t startRunNumber, Int_t endRunNumber, TString ocdbStorage){
244 // make calibration of the chamber status
246 // startRunNumber, endRunNumber - run validity period
247 // ocdbStorage - path to the OCDB storage
248 // - if empty - local storage 'pwd' uesed
249 if (ocdbStorage.Length()<=0) ocdbStorage="local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
252 // 2. extraction of the information
254 if(!AnalyzeChamberStatus()) return;
256 // 3. Append QA plots
258 //MakeDefaultPlots(fVdriftArray,fVdriftArray);
262 // 4. validate OCDB entries
264 if(fSwitchOnValidation==kTRUE && ValidateChamberStatus()==kFALSE) {
265 AliError("TRD Chamber status OCDB parameters not ok!");
272 UpdateOCDBChamberStatus(startRunNumber,endRunNumber,ocdbStorage);
275 //______________________________________________________________________________________________________
276 Bool_t AliTRDPreprocessorOffline::Init(const Char_t* fileName){
278 // read the calibration used during the reconstruction
281 if(ReadVdriftT0Global(fileName)) {
283 TString nameph = fPH2d->GetTitle();
284 fFirstRunVdriftUsed = GetFirstRun(nameph);
285 fVersionVdriftUsed = GetVersion(nameph);
286 fSubVersionVdriftUsed = GetSubVersion(nameph);
288 //printf("Found Version %d, Subversion %d for vdrift\n",fVersionVdriftUsed,fSubVersionVdriftUsed);
292 if(ReadGainGlobal(fileName)) {
294 TString namech = fCH2d->GetTitle();
295 fVersionGainUsed = GetVersion(namech);
296 fSubVersionGainUsed = GetSubVersion(namech);
298 //printf("Found Version %d, Subversion %d for gain\n",fVersionGainUsed,fSubVersionGainUsed);
305 //___________________________________________________________________________________________________________________
307 Bool_t AliTRDPreprocessorOffline::ReadGainGlobal(const Char_t* fileName){
309 // read calibration entries from file
311 if(fCH2d) return kTRUE;
312 TFile fcalib(fileName);
313 TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
315 TH2I *ch2d = (TH2I *) array->FindObject("CH2d");
316 if(!ch2d) return kFALSE;
317 fCH2d = (TH2I*)ch2d->Clone();
318 //fNEvents = (TH1I *) array->FindObject("NEvents");
319 //fAbsoluteGain = (TH2F *) array->FindObject("AbsoluteGain");
321 TH2I *ch2d = (TH2I *) fcalib.Get("CH2d");
322 if(!ch2d) return kFALSE;
323 fCH2d = (TH2I*)ch2d->Clone();
324 //fNEvents = (TH1I *) fcalib.Get("NEvents");
325 //fAbsoluteGain = (TH2F *) fcalib.Get("AbsoluteGain");
327 fCH2d->SetDirectory(0);
328 //printf("title of CH2d %s\n",fCH2d->GetTitle());
333 //_________________________________________________________________________________________________________________
335 Bool_t AliTRDPreprocessorOffline::ReadVdriftT0Global(const Char_t* fileName){
337 // read calibration entries from file
339 if(fPH2d) return kTRUE;
340 TFile fcalib(fileName);
341 TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
343 TProfile2D *ph2d = (TProfile2D *) array->FindObject("PH2d");
344 if(!ph2d) return kFALSE;
345 fPH2d = (TProfile2D*)ph2d->Clone();
346 //fNEvents = (TH1I *) array->FindObject("NEvents");
348 TProfile2D *ph2d = (TProfile2D *) fcalib.Get("PH2d");
349 if(!ph2d) return kFALSE;
350 fPH2d = (TProfile2D*)ph2d->Clone();
351 //fNEvents = (TH1I *) fcalib.Get("NEvents");
353 fPH2d->SetDirectory(0);
354 //printf("title of PH2d %s\n",fPH2d->GetTitle());
359 //___________________________________________________________________________________________________________________
361 Bool_t AliTRDPreprocessorOffline::ReadVdriftLinearFitGlobal(const Char_t* fileName){
363 // read calibration entries from file
365 if(fAliTRDCalibraVdriftLinearFit) return kTRUE;
366 TFile fcalib(fileName);
367 TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
369 fAliTRDCalibraVdriftLinearFit = (AliTRDCalibraVdriftLinearFit *) array->FindObject("AliTRDCalibraVdriftLinearFit");
370 //fNEvents = (TH1I *) array->FindObject("NEvents");
372 fAliTRDCalibraVdriftLinearFit = (AliTRDCalibraVdriftLinearFit *) fcalib.Get("AliTRDCalibraVdriftLinearFit");
373 //fNEvents = (TH1I *) fcalib.Get("NEvents");
375 if(!fAliTRDCalibraVdriftLinearFit) {
376 //printf("No AliTRDCalibraVdriftLinearFit\n");
382 //_____________________________________________________________________________________________________________
384 Bool_t AliTRDPreprocessorOffline::ReadPRFGlobal(const Char_t* fileName){
386 // read calibration entries from file
388 if(fPRF2d) return kTRUE;
389 TFile fcalib(fileName);
390 TObjArray * array = (TObjArray*)fcalib.Get(fNameList);
392 TProfile2D *prf2d = (TProfile2D *) array->FindObject("PRF2d");
393 if(!prf2d) return kFALSE;
394 fPRF2d = (TProfile2D*)prf2d->Clone();
395 //fNEvents = (TH1I *) array->FindObject("NEvents");
397 TProfile2D *prf2d = (TProfile2D *) fcalib.Get("PRF2d");
398 if(!prf2d) return kFALSE;
399 fPRF2d = (TProfile2D*)prf2d->Clone();
400 //fNEvents = (TH1I *) fcalib.Get("NEvents");
402 fPRF2d->SetDirectory(0);
403 //printf("title of PRF2d %s\n",fPRF2d->GetTitle());
408 //__________________________________________________________________________________________________________
410 Bool_t AliTRDPreprocessorOffline::AnalyzeGain(){
412 // Analyze gain - produce the calibration objects
415 AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
416 calibra->SetMinEntries(fMinStatsGain); // If there is less than 1000 entries in the histo: no fit
417 calibra->AnalyseCH(fCH2d);
419 Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(0))
420 + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(0));
421 Int_t nbfit = calibra->GetNumberFit();
422 Int_t nbE = calibra->GetNumberEnt();
426 Bool_t meanother = kFALSE;
429 (nbfit >= 0.5*nbE) && (nbE > 30)) {
430 // create the cal objects
431 if(!fCalDetGainUsed) {
432 calibra->PutMeanValueOtherVectorFit(1,kTRUE);
435 TObjArray object = calibra->GetVectorFit();
436 AliTRDCalDet *calDetGain = calibra->CreateDetObjectGain(&object,meanother);
437 TH1F *coefGain = calDetGain->MakeHisto1DAsFunctionOfDet();
438 // Put them in the array
439 fCalibObjects->AddAt(calDetGain,kGain);
440 fPlots->AddAt(coefGain,kGain);
445 calibra->ResetVectorFit();
450 //_____________________________________________________________________________________________________
451 Bool_t AliTRDPreprocessorOffline::AnalyzeVdriftT0(){
453 // Analyze VdriftT0 - produce the calibration objects
456 AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
457 calibra->SetMinEntries(fMinStatsVdriftT0PH); // If there is less than 1000 entries in the histo: no fit
458 calibra->AnalysePH(fPH2d);
460 Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(1))
461 + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(1));
462 Int_t nbfit = calibra->GetNumberFit();
463 Int_t nbfitSuccess = calibra->GetNumberFitSuccess();
464 Int_t nbE = calibra->GetNumberEnt();
466 //printf("nbtg %d, nbfit %d, nbE %d, nbfitSuccess %d\n",nbtg,nbfit,nbE,nbfitSuccess);
470 (nbfit >= 0.5*nbE) && (nbE > 30) && (nbfitSuccess > 30)) {
471 //printf("Pass the cut for VdriftT0\n");
472 // create the cal objects
473 calibra->RemoveOutliers(1,kFALSE);
474 calibra->PutMeanValueOtherVectorFit(1,kFALSE);
475 calibra->RemoveOutliers2(kFALSE);
476 calibra->PutMeanValueOtherVectorFit2(1,kFALSE);
478 TObjArray object = calibra->GetVectorFit();
479 AliTRDCalDet *calDetVdrift = calibra->CreateDetObjectVdrift(&object);
480 TH1F *coefVdriftPH = calDetVdrift->MakeHisto1DAsFunctionOfDet();
481 AliTRDCalPad *calPadVdrift = (AliTRDCalPad *)calibra->CreatePadObjectVdrift(&object,calDetVdrift);
482 TH1F *coefPadVdrift = calPadVdrift->MakeHisto1D();
483 object = calibra->GetVectorFit2();
484 AliTRDCalDet *calDetT0 = calibra->CreateDetObjectT0(&object);
485 TH1F *coefT0 = calDetT0->MakeHisto1DAsFunctionOfDet();
486 AliTRDCalPad *calPadT0 = (AliTRDCalPad *)calibra->CreatePadObjectT0(&object,calDetT0);
487 TH1F *coefPadT0 = calPadT0->MakeHisto1D();
488 // Put them in the array
489 fCalibObjects->AddAt(calDetT0,kT0PHDet);
490 fCalibObjects->AddAt(calDetVdrift,kVdriftPHDet);
491 fCalibObjects->AddAt(calPadT0,kT0PHPad);
492 fCalibObjects->AddAt(calPadVdrift,kVdriftPHPad);
493 fPlots->AddAt(coefVdriftPH,kVdriftPHDet);
494 fPlots->AddAt(coefT0,kT0PHDet);
495 fPlots->AddAt(coefPadVdrift,kVdriftPHPad);
496 fPlots->AddAt(coefPadT0,kT0PHPad);
500 calibra->ResetVectorFit();
505 //____________________________________________________________________________________________________________________
506 Bool_t AliTRDPreprocessorOffline::AnalyzeVdriftLinearFit(){
508 // Analyze vdrift linear fit - produce the calibration objects
511 //printf("Analyse linear fit\n");
514 AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
515 calibra->SetMinEntries(fMinStatsVdriftLinear); // If there is less than 1000 entries in the histo: no fit
516 //printf("Fill PE Array\n");
517 fAliTRDCalibraVdriftLinearFit->FillPEArray();
518 //printf("AliTRDCalibraFit\n");
519 calibra->AnalyseLinearFitters(fAliTRDCalibraVdriftLinearFit);
523 Int_t nbfit = calibra->GetNumberFit();
524 Int_t nbE = calibra->GetNumberEnt();
529 if ((nbfit >= 0.5*nbE) && (nbE > 30)) {
530 // create the cal objects
531 //calibra->RemoveOutliers(1,kTRUE);
532 calibra->PutMeanValueOtherVectorFit(1,kTRUE);
533 //calibra->RemoveOutliers2(kTRUE);
534 calibra->PutMeanValueOtherVectorFit2(1,kTRUE);
536 TObjArray object = calibra->GetVectorFit();
537 AliTRDCalDet *calDetVdrift = calibra->CreateDetObjectVdrift(&object,kTRUE);
538 TH1F *coefDriftLinear = calDetVdrift->MakeHisto1DAsFunctionOfDet();
539 object = calibra->GetVectorFit2();
540 AliTRDCalDet *calDetLorentz = calibra->CreateDetObjectLorentzAngle(&object);
541 TH1F *coefLorentzAngle = calDetLorentz->MakeHisto1DAsFunctionOfDet();
542 // Put them in the array
543 fCalibObjects->AddAt(calDetVdrift,kVdriftLinear);
544 fCalibObjects->AddAt(calDetLorentz,kLorentzLinear);
545 fPlots->AddAt(coefDriftLinear,kVdriftLinear);
546 fPlots->AddAt(coefLorentzAngle,kLorentzLinear);
551 calibra->ResetVectorFit();
556 //________________________________________________________________________________________________________________
558 Bool_t AliTRDPreprocessorOffline::AnalyzePRF(){
560 // Analyze PRF - produce the calibration objects
563 AliTRDCalibraFit *calibra = AliTRDCalibraFit::Instance();
564 calibra->SetMinEntries(fMinStatsPRF); // If there is less than 1000 entries in the histo: no fit
565 calibra->AnalysePRFMarianFit(fPRF2d);
567 Int_t nbtg = 6*4*18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb0(2))
568 + 6* 18*((Int_t) ((AliTRDCalibraMode *)calibra->GetCalibraMode())->GetDetChamb2(2));
569 Int_t nbfit = calibra->GetNumberFit();
570 Int_t nbE = calibra->GetNumberEnt();
576 (nbfit >= 0.95*nbE) && (nbE > 30)) {
577 // create the cal objects
578 TObjArray object = calibra->GetVectorFit();
579 AliTRDCalPad *calPadPRF = (AliTRDCalPad*) calibra->CreatePadObjectPRF(&object);
580 TH1F *coefPRF = calPadPRF->MakeHisto1D();
581 // Put them in the array
582 fCalibObjects->AddAt(calPadPRF,kPRF);
583 fPlots->AddAt(coefPRF,kPRF);
588 calibra->ResetVectorFit();
594 //_____________________________________________________________________________
595 Bool_t AliTRDPreprocessorOffline::AnalyzeChamberStatus()
598 // Produce AliTRDCalChamberStatus out of calibration results
601 // set up AliTRDCalChamberStatus
602 AliTRDCalChamberStatus *CalChamberStatus = new AliTRDCalChamberStatus();
603 for(Int_t det = 0; det < 540; det++) CalChamberStatus->SetStatus(det,1);
605 // get calibration objects
606 AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain);
607 AliTRDCalDet *calDetVDrift = (AliTRDCalDet *) fCalibObjects->At(kVdriftLinear);
610 if((!calDetGain) || (!calDetVDrift) || (!fCH2d)) return kFALSE;
613 Double_t gainmean = calDetGain->GetMean();
614 Double_t vdriftmean = calDetVDrift->GetMean();
616 Double_t gainrms = calDetGain->GetRMSRobust();
617 Double_t vdriftrms = calDetVDrift->GetRMSRobust();
619 //printf("Gain mean: %f, rms: %f\n",gainmean,gainrms);
620 //printf("Vdrift mean: %f, rms: %f\n",vdriftmean,vdriftrms);
623 if((TMath::Abs(gainrms) < 0.001) || (TMath::Abs(vdriftrms) < 0.001)) return kFALSE;
625 // mask chambers with empty gain entries
627 for (Int_t idet = 0; idet < 540; idet++) {
630 TH1I *projch = (TH1I *) fCH2d->ProjectionX("projch",idet+1,idet+1,(Option_t *)"e");
631 Double_t entries = projch->GetEntries();
634 Double_t gain = calDetGain->GetValue(idet);
637 Double_t vdrift = calDetVDrift->GetValue(idet);
641 TMath::Abs(gainmean-gain) > (15.0*gainrms) ||
642 TMath::Abs(vdriftmean-vdrift) > (15.0*vdriftrms)) {
644 //printf(" chamber det %03d masked \n",idet);
645 //printf(" gainmean %f and gain %f, gainrms %f \n",gainmean,gain,gainrms);
646 //printf(" vdriftmean %f and vdrift %f, vdriftrms %f \n",vdriftmean,vdrift,vdriftrms);
647 CalChamberStatus->SetStatus(idet,AliTRDCalChamberStatus::kMasked);
652 // installed supermodules+1 -> abort
653 if(counter > (7+1)*30) {
654 printf("ERROR: more than one SM to be masked!! \n Abort...\n");
655 if(projch) delete projch;
665 for(Int_t sm=0; sm < 18; sm++) {
667 for(Int_t det = 0; det < 30; det++){
668 Int_t detector = sm*30+det;
669 if(CalChamberStatus->IsMasked(detector)) counter++;
672 for(Int_t det = 0; det < 30; det++){
673 Int_t detector = sm*30+det;
674 CalChamberStatus->SetStatus(detector,AliTRDCalChamberStatus::kInstalled);
679 fCalibObjects->AddAt(CalChamberStatus,kChamberStatus);
685 //________________________________________________________________________________________________
686 void AliTRDPreprocessorOffline::CorrectFromDetGainUsed() {
688 // Correct from the gas gain used afterwards
690 AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain);
691 if(!calDetGain) return;
697 for(Int_t det = 0; det < 540; det++) {
699 Float_t gaininit = fCalDetGainUsed->GetValue(det);
700 Float_t gainout = calDetGain->GetValue(det);
703 if(TMath::Abs(gainout-1.0) > 0.000001) {
704 mean += (gaininit*gainout);
708 if(nbdet > 0) mean = mean/nbdet;
710 for(Int_t det = 0; det < 540; det++) {
712 Float_t gaininit = fCalDetGainUsed->GetValue(det);
713 Float_t gainout = calDetGain->GetValue(det);
715 if(TMath::Abs(gainout-1.0) > 0.000001) {
716 Double_t newgain = gaininit*gainout;
717 if(newgain < 0.1) newgain = 0.1;
718 if(newgain > 1.9) newgain = 1.9;
719 calDetGain->SetValue(det,newgain);
722 Double_t newgain = mean;
723 if(newgain < 0.1) newgain = 0.1;
724 if(newgain > 1.9) newgain = 1.9;
725 calDetGain->SetValue(det,newgain);
731 //________________________________________________________________________________________________
732 void AliTRDPreprocessorOffline::CorrectFromDetVdriftUsed() {
734 // Correct from the drift velocity
737 //printf("Correct for vdrift\n");
739 AliTRDCalDet *calDetGain = (AliTRDCalDet *) fCalibObjects->At(kGain);
740 if(!calDetGain) return;
742 Int_t detVdrift = kVdriftPHDet;
743 if(fMethodSecond) detVdrift = kVdriftLinear;
744 AliTRDCalDet *calDetVdrift = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
745 if(!calDetVdrift) return;
749 for(Int_t det = 0; det < 540; det++) {
751 Float_t vdriftinit = fCalDetVdriftUsed->GetValue(det);
752 Float_t vdriftout = calDetVdrift->GetValue(det);
754 Float_t gain = calDetGain->GetValue(det);
755 if(vdriftout > 0.0) gain = gain*vdriftinit/vdriftout;
756 if(gain < 0.1) gain = 0.1;
757 if(gain > 1.9) gain = 1.9;
758 calDetGain->SetValue(det,gain);
764 //_________________________________________________________________________________________________________________
765 void AliTRDPreprocessorOffline::UpdateOCDBGain(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
770 AliCDBMetaData *metaData= new AliCDBMetaData();
771 metaData->SetObjectClassName("AliTRDCalDet");
772 metaData->SetResponsible("Raphaelle Bailhache");
773 metaData->SetBeamPeriod(1);
775 AliCDBId id1("TRD/Calib/ChamberGainFactor", startRunNumber, endRunNumber);
776 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
777 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain);
778 if(calDet) gStorage->Put(calDet, id1, metaData);
782 //___________________________________________________________________________________________________________________
783 void AliTRDPreprocessorOffline::UpdateOCDBVdrift(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
788 Int_t detVdrift = kVdriftPHDet;
790 if(fMethodSecond) detVdrift = kVdriftLinear;
792 AliCDBMetaData *metaData= new AliCDBMetaData();
793 metaData->SetObjectClassName("AliTRDCalDet");
794 metaData->SetResponsible("Raphaelle Bailhache");
795 metaData->SetBeamPeriod(1);
797 AliCDBId id1("TRD/Calib/ChamberVdrift", startRunNumber, endRunNumber);
798 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
799 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
800 if(calDet) gStorage->Put(calDet, id1, metaData);
806 AliCDBMetaData *metaDataPad= new AliCDBMetaData();
807 metaDataPad->SetObjectClassName("AliTRDCalPad");
808 metaDataPad->SetResponsible("Raphaelle Bailhache");
809 metaDataPad->SetBeamPeriod(1);
811 AliCDBId id1Pad("TRD/Calib/LocalVdrift", startRunNumber, endRunNumber);
812 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad);
813 if(calPad) gStorage->Put(calPad, id1Pad, metaDataPad);
818 //________________________________________________________________________________________________________________________
819 void AliTRDPreprocessorOffline::UpdateOCDBT0(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
824 AliCDBMetaData *metaData= new AliCDBMetaData();
825 metaData->SetObjectClassName("AliTRDCalDet");
826 metaData->SetResponsible("Raphaelle Bailhache");
827 metaData->SetBeamPeriod(1);
829 AliCDBId id1("TRD/Calib/ChamberT0", startRunNumber, endRunNumber);
830 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
831 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet);
832 if(calDet) gStorage->Put(calDet, id1, metaData);
836 AliCDBMetaData *metaDataPad= new AliCDBMetaData();
837 metaDataPad->SetObjectClassName("AliTRDCalPad");
838 metaDataPad->SetResponsible("Raphaelle Bailhache");
839 metaDataPad->SetBeamPeriod(1);
841 AliCDBId id1Pad("TRD/Calib/LocalT0", startRunNumber, endRunNumber);
842 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad);
843 if(calPad) gStorage->Put(calPad, id1Pad, metaDataPad);
848 //_________________________________________________________________________________________________________________
849 void AliTRDPreprocessorOffline::UpdateOCDBPRF(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
854 AliCDBMetaData *metaData= new AliCDBMetaData();
855 metaData->SetObjectClassName("AliTRDCalPad");
856 metaData->SetResponsible("Raphaelle Bailhache");
857 metaData->SetBeamPeriod(1);
859 AliCDBId id1("TRD/Calib/PRFWidth", startRunNumber, endRunNumber);
860 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
861 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF);
862 if(calPad) gStorage->Put(calPad, id1, metaData);
866 //_________________________________________________________________________________________________________________
867 void AliTRDPreprocessorOffline::UpdateOCDBChamberStatus(Int_t startRunNumber, Int_t endRunNumber, const Char_t *storagePath){
872 AliCDBMetaData *metaData= new AliCDBMetaData();
873 metaData->SetObjectClassName("AliTRDCalChamberStatus");
874 metaData->SetResponsible("Raphaelle Bailhache");
875 metaData->SetBeamPeriod(1);
877 AliCDBId id1("TRD/Calib/ChamberStatus", startRunNumber, endRunNumber);
878 AliCDBStorage * gStorage = AliCDBManager::Instance()->GetStorage(storagePath);
879 AliTRDCalChamberStatus *calChamberStatus = (AliTRDCalChamberStatus *) fCalibObjects->At(kChamberStatus);
880 if(calChamberStatus) gStorage->Put(calChamberStatus, id1, metaData);
884 //__________________________________________________________________________________________________________________________
885 Bool_t AliTRDPreprocessorOffline::ValidateGain() const {
887 // Validate OCDB entry
890 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kGain);
892 Double_t mean = calDet->GetMean();
893 Double_t rms = calDet->GetRMSRobust();
894 if((mean > 0.2) && (mean < 1.4) && (rms < 0.5)) return kTRUE;
895 //if((mean > 0.2) && (mean < 1.4)) return kTRUE;
902 //__________________________________________________________________________________________________________________________
903 Bool_t AliTRDPreprocessorOffline::ValidateVdrift(){
908 Int_t detVdrift = kVdriftPHDet;
911 if(fMethodSecond) detVdrift = kVdriftLinear;
913 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(detVdrift);
915 Double_t mean = calDet->GetMean();
916 Double_t rms = calDet->GetRMSRobust();
917 //printf("Vdrift::mean %f, rms %f\n",mean,rms);
918 if(!((mean > 1.0) && (mean < 2.0) && (rms < 0.5))) ok = kFALSE;
923 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kVdriftPHPad);
925 Double_t mean = calPad->GetMean();
926 Double_t rms = calPad->GetRMS();
927 //printf("Vdrift::meanpad %f, rmspad %f\n",mean,rms);
928 if(!((mean > 0.9) && (mean < 1.1) && (rms < 0.6))) ok = kFALSE;
936 //__________________________________________________________________________________________________________________________
937 Bool_t AliTRDPreprocessorOffline::ValidateT0(){
942 AliTRDCalDet *calDet = (AliTRDCalDet *) fCalibObjects->At(kT0PHDet);
943 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kT0PHPad);
944 if(calDet && calPad) {
945 Double_t meandet = calDet->GetMean();
946 Double_t rmsdet = calDet->GetRMSRobust();
947 Double_t meanpad = calPad->GetMean();
948 //Double_t rmspad = calPad->GetRMS();
949 //printf("T0::minimum %f, rmsdet %f,meanpad %f, rmspad %f\n",meandet,rmsdet,meanpad,rmspad);
950 if((meandet > -1.5) && (meandet < 5.0) && (rmsdet < 4.0) && (meanpad < 5.0) && (meanpad > -0.5)) return kTRUE;
956 //__________________________________________________________________________________________________________________________
957 Bool_t AliTRDPreprocessorOffline::ValidatePRF() const{
962 AliTRDCalPad *calPad = (AliTRDCalPad *) fCalibObjects->At(kPRF);
964 Double_t meanpad = calPad->GetMean();
965 Double_t rmspad = calPad->GetRMS();
966 //printf("PRF::meanpad %f, rmspad %f\n",meanpad,rmspad);
967 if((meanpad < 1.0) && (rmspad < 0.8)) return kTRUE;
974 //__________________________________________________________________________________________________________________________
975 Bool_t AliTRDPreprocessorOffline::ValidateChamberStatus() const{
980 AliTRDCalChamberStatus *calChamberStatus = (AliTRDCalChamberStatus *) fCalibObjects->At(kChamberStatus);
981 if(calChamberStatus) {
982 Int_t detectormasked = 0;
983 for(Int_t det = 0; det < 540; det++) {
984 if(calChamberStatus->IsMasked(det)) detectormasked++;
986 //printf("Number of chambers masked %d\n",detectormasked);
987 if(detectormasked > 29) return kFALSE;
993 //_____________________________________________________________________________
994 Int_t AliTRDPreprocessorOffline::GetVersion(TString name) const
997 // Get version from the title
1001 const Char_t *version = "Ver";
1002 if(!strstr(name.Data(),version)) return -1;
1003 const Char_t *after = "Subver";
1004 if(!strstr(name.Data(),after)) return -1;
1006 for(Int_t ver = 0; ver < 999999999; ver++) {
1008 TString vertry(version);
1012 //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
1014 if(strstr(name.Data(),vertry.Data())) return ver;
1022 //_____________________________________________________________________________
1023 Int_t AliTRDPreprocessorOffline::GetSubVersion(TString name) const
1026 // Get subversion from the title
1030 const Char_t *subversion = "Subver";
1031 if(!strstr(name.Data(),subversion)) return -1;
1032 const Char_t *after = "FirstRun";
1033 if(!strstr(name.Data(),after)) {
1036 if(!strstr(name.Data(),after)) return -1;
1039 for(Int_t ver = 0; ver < 999999999; ver++) {
1041 TString vertry(subversion);
1045 //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
1047 if(strstr(name.Data(),vertry.Data())) return ver;
1055 //_____________________________________________________________________________
1056 Int_t AliTRDPreprocessorOffline::GetFirstRun(TString name) const
1059 // Get first run from the title
1063 const Char_t *firstrun = "FirstRun";
1064 if(!strstr(name.Data(),firstrun)) return -1;
1065 const Char_t *after = "Nz";
1066 if(!strstr(name.Data(),after)) return -1;
1069 for(Int_t ver = 0; ver < 999999999; ver++) {
1071 TString vertry(firstrun);
1075 //printf("vertry %s and name %s\n",vertry.Data(),name.Data());
1077 if(strstr(name.Data(),vertry.Data())) return ver;