]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/macrosSDD/MakeOCDBCorrectionFromOutputADCCalib.C
commented define _ClusterTopology_ - to be used only for the special productions
[u/mrichter/AliRoot.git] / ITS / macrosSDD / MakeOCDBCorrectionFromOutputADCCalib.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TH1F.h>
3 #include <TH2F.h>
4 #include <TF1.h>
5 #include <TPad.h>
6 #include <TGraphErrors.h>
7 #include <TROOT.h>
8 #include <TFile.h>
9 #include <TTree.h>
10 #include <TGrid.h>
11 #include <TGridResult.h>
12 #include <TMath.h>
13 #include <TCanvas.h>
14 #include <TStyle.h>
15 #include <TLatex.h>
16 #include <AliCDBEntry.h>
17 #include <AliITSresponseSDD.h>
18
19 #endif
20
21 // Inputs: Output of MakeSDDCalib.C + AliITSresponseSDD Object used to produce the cpass1 (default Run192772_192779_v0_s0.root)
22 // Outputs: new ADCtoKeV and ADCvsDriftTimeHistos
23 // Launch as: aliroot MakeCorrectionFromOutputCalib.C+
24
25 void MakeOCDBCorrectionFromOutputADCCalib(TString fname=".",TString OCDBname="Run192772_192779_v0_s0")
26 {
27   Printf("Opening file CalibResults.root in %s/",fname.Data());
28   TFile *fin = new TFile(Form("%s/SDDADCCalibResults.root",fname.Data()));
29     
30   TH1F *hmpvModpar0=(TH1F*)fin->Get("hmpvModpar0");
31   TH1F *hmpvModpar1=(TH1F*)fin->Get("hmpvModpar1");
32   Float_t s=0.0101;//Default ADC vs drift time
33   // getting ADC2keV
34   TFile* fr=TFile::Open(Form("%s.root",OCDBname.Data()));
35   AliCDBEntry* e=(AliCDBEntry*)fr->Get("AliCDBEntry");
36   AliITSresponseSDD* r=(AliITSresponseSDD*)e->GetObject();
37   TH1F* hADCtokeV=new TH1F("hADCtokeV","",260,239.5,499.5);
38   TH1F* hADCvsDriftTime=new TH1F("hADCvsDriftTime","",260,239.5,499.5);
39   r->ls();
40     
41   for(Int_t iMod=240; iMod<500; iMod++){
42     Float_t ak=r->GetADCtokeV(iMod);
43     Float_t adcdrtime=r->GetADCvsDriftTime(iMod);
44     printf("mod %d\nadc->keV=%f ",iMod,ak);
45     hADCtokeV->SetBinContent(iMod-240+1,ak);
46     hADCvsDriftTime->SetBinContent(iMod-240+1,adcdrtime);
47     //calculate new ADC2KeV
48     Double_t kprim=hmpvModpar0->GetBinContent(hmpvModpar0->FindBin(iMod));
49     if(kprim==0)kprim=84;
50     if(iMod==376)kprim=84;
51     Double_t Corr=kprim*ak/84;
52     printf("mpv=%f | Corr(adc->keV*mpv/84)=%f\n",kprim,Corr);
53     hmpvModpar0->SetBinContent(hmpvModpar0->FindBin(iMod),Corr);
54     //Calculate new ADCvsDrTime
55     Double_t sprim=hmpvModpar1->GetBinContent(hmpvModpar1->FindBin(iMod));
56     if(iMod==376){
57       sprim=0;
58     }
59     Double_t Corr2=adcdrtime-(sprim*ak/2);
60     hmpvModpar1->SetBinContent(hmpvModpar1->FindBin(iMod),Corr2);
61     printf("ADCvsDriftTime=%f | Corr(ADCvsDriftTime_old-(ADCvsDriftTime_new*adc->keV/2))=%f\n",adcdrtime,Corr2);
62   }
63     
64   hmpvModpar0->GetYaxis()->SetTitle("K^{*}=KK'/84");
65   hmpvModpar1->GetYaxis()->SetTitle("S^{*}=S-S'K/2");
66   hmpvModpar0->SetTitle("hADCtokeV - New");
67   hmpvModpar0->SetName("hADCtokeV");
68   hADCtokeV->SetTitle("hADCtokeV - Run166530_999999999_v4_s0");
69   hADCtokeV->SetLineColor(2);
70   TCanvas *c0=new TCanvas("c0","c0");
71   c0->cd(1);
72   hmpvModpar0->DrawClone("HIST");
73   hADCtokeV->DrawClone("HISTSAME");
74   gPad->BuildLegend();
75     
76     
77   hmpvModpar1->SetTitle("hADCvsDriftTime - New");
78   hmpvModpar1->SetName("hADCvsDriftTime");
79   hADCvsDriftTime->SetTitle("hADCvsDriftTime - Run166530_999999999_v4_s0");
80   hADCvsDriftTime->SetLineColor(2);
81   TCanvas *c1=new TCanvas("c1","c1");
82   hmpvModpar1->DrawClone("HIST");
83   hADCvsDriftTime->DrawClone("HISTSAME");
84   c1->BuildLegend();
85   TFile *out=new TFile(Form("%s/CorrectiondEdxSDD_%s_%s.root",fname.Data(),fname.Data(),OCDBname.Data()),"recreate");
86   hmpvModpar0->Write();
87   hmpvModpar1->Write();
88   out->Close();
89   delete out;
90     
91 } //end main
92