]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/Upgrade/macros/makeResidualSpaceChargeOCDB.C
f6c5eddbf1831732fc530dc71198c103ca3f15d8
[u/mrichter/AliRoot.git] / TPC / Upgrade / macros / makeResidualSpaceChargeOCDB.C
1 //
2 // macro to create a residual OCDB 
3 //
4
5 /*
6 //
7 //
8 //
9 .x $ALICE_ROOT/TPC/Upgrade/macros/ConfigOCDB.C(1)
10 .L $ALICE_ROOT/TPC/Upgrade/macros/makeResidualSpaceChargeOCDB.C
11
12 Example usage:
13
14 ln -sf /hera/alice/wiechula/Upgrade/LUTs_fluctuation/dir1/SpaceChargeFluc10_1.lookup.root current.root
15 ln -sf /hera/alice/wiechula/Upgrade/LUTs_fluctuation/average/SpaceChargeFluc0_1.lookup.root  mean.root
16 ln -sf $ALICE_ROOT/OCDB/TPC/Calib/Correction/Run0_999999999_v0_s2.root ocdb.root
17 aliroot -b -q $ALICE_ROOT/TPC/Upgrade/macros/makeResidualSpaceChargeOCDB.C
18
19 */
20
21 void makeResidualSpaceChargeOCDBLookup(const char *ocdbInName="ocdb.root",const char *scCurrentName="current.root", const char *scMeanName="mean.root"){
22   //
23   // Macro to create a clone original  OCDB entry with space point distortion and add there space point
24   // distortions cause by resifual sapce charge
25   // Output is stored by default in the current directory.
26   // 
27
28   // Parameters to Specify:
29   //    ocdbInName      : path to the original OCDB entry
30   //    scCurrentName   : path to the fluctuetaed distortion map
31   //    scMeanName      : path to the mean distortion map
32   //    
33   
34   /*
35     scCurrentName="/hera/alice/wiechula/Upgrade/LUTs_fluctuation/dir1/SpaceChargeFluc10_1.lookup.root";
36     scMeanName="/hera/alice/wiechula/Upgrade/LUTs_fluctuation/average/SpaceChargeFluc0_1.lookup.root";
37     ocdbInName="$ALICE_ROOT/OCDB/TPC/Calib/Correction/Run0_999999999_v0_s2.root"
38   */
39   TFile * finCurrent = TFile::Open(scCurrentName);
40   TFile * finMean = TFile::Open(scMeanName);
41   //
42   AliTPCCorrectionLookupTable *spaceChargeCurrent=  (AliTPCCorrectionLookupTable *)finCurrent->Get("map");
43   AliTPCCorrectionLookupTable *spaceChargeMean   =  (AliTPCCorrectionLookupTable *)finMean->Get("map");
44   AliTPCInverseCorrection * spaceChargeInverseMean = new AliTPCInverseCorrection(spaceChargeMean);
45   //
46   TObjArray * arraySC = new TObjArray(2);
47   arraySC->AddAt(spaceChargeCurrent,0);
48   arraySC->AddAt(spaceChargeMean,1);
49   AliTPCComposedCorrection *corrSC = new AliTPCComposedCorrection(arraySC,AliTPCComposedCorrection::kParallel);
50   AliTPCComposedCorrection *corrSCW = new AliTPCComposedCorrection(arraySC,AliTPCComposedCorrection::kParallel);
51   AliTPCComposedCorrection::AddVisualCorrection(spaceChargeCurrent,1);
52   AliTPCComposedCorrection::AddVisualCorrection(spaceChargeMean,2);
53   AliTPCComposedCorrection::AddVisualCorrection(spaceChargeInverseMean,3);
54   AliTPCComposedCorrection::AddVisualCorrection(corrSC,4);
55   AliTPCComposedCorrection::AddVisualCorrection(corrSCW,5);
56   //
57   //
58   TLinearFitter fitterR(2,"pol1");
59   for (Int_t ipoint =0; ipoint<10000; ipoint++){
60     Double_t z0 = -250+gRandom->Rndm()*500;
61     Double_t r0   = 85+gRandom->Rndm()*(245-85.);
62     Double_t phi0 = gRandom->Rndm()*TMath::TwoPi();
63     // some problematic parts to be skipped  - investigated later
64     Double_t xvalue[2]={0};
65     if (TMath::Abs(AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,0,1))>50) continue;
66     if (TMath::Abs(AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,0,2))>50) continue;
67     if (TMath::Abs(AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,1,1))>20) continue;
68     if (TMath::Abs(AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,1,2))>20) continue;
69     if (TMath::Abs(AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,2,1))>50) continue;
70     if (TMath::Abs(AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,2,2))>50) continue;
71     //
72     //
73     xvalue[0]=AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,0,2);
74     fitterR.AddPoint(xvalue,AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,0,1));
75   }
76   fitterR->Eval();
77   TVectorD weights(2);
78   weights[0]=1;
79   weights[1]=-1;
80   corrSC->SetWeights( &weights);
81   weights[1]=-fitterR->GetParameter(1);
82   corrSCW->SetWeights( &weights);
83   //
84   //
85   //
86   TF1 *fSC = new TF1("fSC","AliTPCCorrection::GetCorrXYZ(x,0,10,0,4)",85,245);  
87   TF1 *fSCW = new TF1("fSCW","AliTPCCorrection::GetCorrXYZ(x,0,10,0,5)",85,245);
88   fSC->SetLineColor(2);
89   fSC->SetLineColor(4);
90   fSC->Draw();
91   fSCW->Draw("same");
92   gPad->SaveAs("residualMap.pdf");
93   
94   //
95   //
96   //
97   TFile *fileOCDBIn=TFile::Open(ocdbInName);
98   AliCDBEntry   *entry =  ( AliCDBEntry   *)fileOCDBIn->Get("AliCDBEntry");
99   TObjArray * array = (TObjArray*)entry->GetObject();
100   for (Int_t i=0;  i<3;  i++){
101     AliTPCComposedCorrection *corr = ( AliTPCComposedCorrection*)array->At(i);
102     if (corr){
103       TObjArray* corrArray = corr->GetCorrections();
104       corrArray->AddLast(corrSCW);
105     }    
106   }
107   AliCDBMetaData *metaData= new AliCDBMetaData();
108   metaData->SetObjectClassName("TObjArray");
109   metaData->SetResponsible("Marian Ivanov");
110   metaData->SetBeamPeriod(1);
111   metaData->SetAliRootVersion("05-25-01"); //root version
112   metaData->SetComment("Standard+fluctuation");
113   AliCDBId* id1=NULL;
114   id1=new AliCDBId("TPC/Calib/Correction", 0, 9999999999999999);
115   //
116   TString localStorage = "local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
117   pocdbStorage = AliCDBManager::Instance()->GetStorage(localStorage.Data());
118   pocdbStorage->Put(array, (*id1), metaData);
119
120
121 }
122
123
124
125
126
127 void makeResidualSpaceChargeOCDB(const char *ocdbInName="ocdb.root",const char *scCurrentName="current.root", const char *scMeanName="mean.root"){
128   //
129   // Macro to create a clone original  OCDB entry with space point distortion and add there space point
130   // distortions cause by resifual sapce charge
131   // Output is stored by default in the current directory.
132   // 
133
134   // Parameters to Specify:
135   //    ocdbInName      : path to the original OCDB entry
136   //    scCurrentName   : path to the fluctuetaed distortion map
137   //    scMeanName      : path to the mean distortion map
138   //    
139   
140   /*
141     scCurrentName="/hera/alice/miranov/SpaceCharge/Fluctuations/PbPbWithGain/dirmergeAll/dEpsilon20/dir0/SpaceChargeFluc10_1.root";
142     scMeanName="/hera/alice/miranov/SpaceCharge/Fluctuations/PbPbWithGain/dirmergeAll/dEpsilon20/dir0//SpaceChargeFluc0_1.root";
143     ocdbInName="$ALICE_ROOT/OCDB/TPC/Calib/Correction/Run0_999999999_v0_s2.root"
144   */
145   TFile * finCurrent = TFile::Open(scCurrentName);
146   TFile * finMean = TFile::Open(scMeanName);
147   //
148   AliTPCCorrectionLookupTable *spaceChargeCurrent=  (AliTPCCorrectionLookupTable *)finCurrent->Get("DistRef");
149   AliTPCCorrectionLookupTable *spaceChargeMean   =  (AliTPCCorrectionLookupTable *)finMean->Get("DistRef");
150   AliTPCInverseCorrection * spaceChargeInverseMean = new AliTPCInverseCorrection(spaceChargeMean);
151   //
152   TObjArray * arraySC = new TObjArray(2);
153   arraySC->AddAt(spaceChargeCurrent,0);
154   arraySC->AddAt(spaceChargeMean,1);
155   AliTPCComposedCorrection *corrSC = new AliTPCComposedCorrection(arraySC,AliTPCComposedCorrection::kParallel);
156   AliTPCComposedCorrection *corrSCW = new AliTPCComposedCorrection(arraySC,AliTPCComposedCorrection::kParallel);
157   AliTPCComposedCorrection::AddVisualCorrection(spaceChargeCurrent,1);
158   AliTPCComposedCorrection::AddVisualCorrection(spaceChargeMean,2);
159   AliTPCComposedCorrection::AddVisualCorrection(spaceChargeInverseMean,3);
160   AliTPCComposedCorrection::AddVisualCorrection(corrSC,4);
161   AliTPCComposedCorrection::AddVisualCorrection(corrSCW,5);
162   //
163   //
164   TLinearFitter fitterR(2,"pol1");
165   for (Int_t ipoint =0; ipoint<10000; ipoint++){
166     Double_t z0 = -250+gRandom->Rndm()*500;
167     Double_t r0   = 85+gRandom->Rndm()*(245-85.);
168     Double_t phi0 = gRandom->Rndm()*TMath::TwoPi();
169     // some problematic parts to be skipped  - investigated later
170     Double_t xvalue[2]={0};
171     if (TMath::Abs(AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,0,1))>50) continue;
172     if (TMath::Abs(AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,0,2))>50) continue;
173     if (TMath::Abs(AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,1,1))>20) continue;
174     if (TMath::Abs(AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,1,2))>20) continue;
175     if (TMath::Abs(AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,2,1))>50) continue;
176     if (TMath::Abs(AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,2,2))>50) continue;
177     //
178     //
179     xvalue[0]=AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,1,2);
180     fitterR.AddPoint(xvalue,AliTPCCorrection::GetCorrXYZ(r0*TMath::Cos(phi0), r0*TMath::Sin(phi0), z0,1,1));
181   }
182   fitterR->Eval();
183   TVectorD weights(2);
184   weights[0]=1;
185   weights[1]=-1;
186   corrSC->SetWeights( &weights);
187   weights[1]=-fitterR->GetParameter(1);
188   corrSCW->SetWeights( &weights);
189   //
190   //
191   //
192   TF1 *fSC = new TF1("fSC","AliTPCCorrection::GetCorrXYZ(x,0,10,1,4)",85,245);  
193   TF1 *fSCW = new TF1("fSCW","AliTPCCorrection::GetCorrXYZ(x,0,10,1,5)",85,245);
194   fSC->SetLineColor(2);
195   fSC->SetLineColor(4);
196   fSC->GetXaxis()->SetTitle(" R (cm)");
197   fSC->GetYaxis()->SetTitle(" #Delta_{R#phi} (cm)");
198   TCanvas * canvasDist = new TCanvas("canvasDist","canvasDist",600,500);
199   TF1 *f1sc = new TF1("f1sc","[0]+[1]*x+[2]*x**2",85,245);
200   TF1 *f1scw = new TF1("f1scw","[0]+[1]*x+[2]*x**2",85,245);
201
202   {
203     fSC->Draw();
204     fSCW->Draw();
205     fSC->SetMinimum(-1);
206     fSC->SetMaximum(2.);
207     fSC->GetHistogram()->SetLineColor(2);
208     fSCW->GetHistogram()->SetLineColor(4);
209     //
210     f1sc->SetLineColor(2);
211     f1scw->SetLineColor(4);
212     f1sc->SetLineWidth(0.5);
213     f1scw->SetLineWidth(0.5);
214     //
215     fSC->GetHistogram()->Fit(f1sc);
216     fSCW->GetHistogram()->Fit(f1scw);
217     //
218     fSC->GetHistogram()->Draw("");
219     fSCW->GetHistogram()->Draw("same");
220     //
221     TF1 *f1scp=new TF1(*f1sc);
222     TF1 *f1scm=new TF1(*f1sc);
223     f1scp->SetLineStyle(2); f1scm->SetLineStyle(2);
224     f1scp->SetParameter(0, f1sc->GetParameter(0)+ 0.3);
225     f1scm->SetParameter(0, f1sc->GetParameter(0)- 0.3);
226     f1scp->Draw("same");
227     f1scm->Draw("same");
228     //
229     TF1 *f1scwp=new TF1(*f1scw);
230     TF1 *f1scwm=new TF1(*f1scw);
231     f1scwp->SetLineStyle(2); f1scwm->SetLineStyle(2);
232     f1scwp->SetParameter(0, f1scw->GetParameter(0)+ 0.3);
233     f1scwm->SetParameter(0, f1scw->GetParameter(0)- 0.3);
234     f1scwp->Draw("same");
235     f1scwm->Draw("same");
236   }
237   TLegend * legend = new TLegend(0.3,0.5,0.89,0.89,"Residual #Delta_{R#phi} distortions and helix fit");
238   legend->SetBorderSize(0);
239   legend->AddEntry(fSC->GetHistogram(),"Stage 0: #Delta_{R#phic}-#Delta_{R#phimean}");
240   legend->AddEntry(fSCW->GetHistogram(),"Stage 1: #Delta_{R#phic}-k_{fit}#Delta_{R#phimean}");
241   legend->Draw();
242   gPad->SaveAs("residualMapTrackFit.pdf");
243   
244   //
245   //
246   //
247   TFile *fileOCDBIn=TFile::Open(ocdbInName);
248   AliCDBEntry   *entry =  ( AliCDBEntry   *)fileOCDBIn->Get("AliCDBEntry");
249   TObjArray * array = (TObjArray*)entry->GetObject();
250   for (Int_t i=0;  i<3;  i++){
251     AliTPCComposedCorrection *corr = ( AliTPCComposedCorrection*)array->At(i);
252     if (corr){
253       TObjArray* corrArray = corr->GetCorrections();
254       corrArray->AddLast(corrSCW);
255     }    
256   }
257   AliCDBMetaData *metaData= new AliCDBMetaData();
258   metaData->SetObjectClassName("TObjArray");
259   metaData->SetResponsible("Marian Ivanov");
260   metaData->SetBeamPeriod(1);
261   metaData->SetAliRootVersion("05-25-01"); //root version
262   metaData->SetComment("Standard+fluctuation");
263   AliCDBId* id1=NULL;
264   id1=new AliCDBId("TPC/Calib/Correction", 0, 9999999999999999);
265   //
266   TString localStorage = "local://"+gSystem->GetFromPipe("pwd")+"/OCDB";
267   pocdbStorage = AliCDBManager::Instance()->GetStorage(localStorage.Data());
268   pocdbStorage->Put(array, (*id1), metaData);
269
270
271 }
272
273