]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/CalibMacros/CalibAlign.C
Update master to aliroot
[u/mrichter/AliRoot.git] / TPC / CalibMacros / CalibAlign.C
1 /// \file CalibAlign.C
2 ///
3 /// ~~~{.cpp}
4 /// .x ~/UliStyle.C
5 /// .x ~/rootlogon.C
6 /// gSystem->Load("libSTAT");
7 /// gSystem->Load("libANALYSIS");
8 /// gSystem->Load("libTPCcalib"); 
9 /// gSystem->Load("libSTAT");
10 /// gSystem->AddIncludePath("-I$ALICE_ROOT/TPC/macros");
11 /// 
12 ///
13 ///
14 /// gROOT->LoadMacro("$ALICE_ROOT/TPC/macros/AliXRDPROOFtoolkit.cxx+")
15 /// AliXRDPROOFtoolkit tool; 
16 /// TChain * chainCosmic = tool.MakeChain("cosmic.txt","Track0",0,1000000);
17 /// chainCosmic->Lookup();
18 /// chainCosmic->SetProof(kTRUE);
19 /// TChain * chainTr = tool.MakeChain("align.txt","Tracklet",0,10200);
20 /// chainTr->Lookup();
21 /// //chainTr->SetProof(kTRUE);
22 ///
23 /// .L $ALICE_ROOT/TPC/CalibMacros/CalibAlign.C
24 /// SetAlias();
25 /// InitCutsAlign();
26 /// MakeAlign();
27 /// ~~~
28   
29 /*
30 #include "TMath.h"
31 #include "TFile.h"
32 #include "TLinearFitter.h"
33 #include "TChain.h"
34 #include "TTreeStream.h"
35 #include "TStatToolkit.h"
36 #include "TH2F.h"
37 #include "TLegend.h"
38 #include "TCut.h"
39 #include "TPad.h"
40 #include "AliTPCcalibAlign.h"
41 */
42
43 AliTPCcalibAlign align;
44 TChain * chainTr;
45
46 void SetAlias(){
47   ///
48
49   chainTr->SetAlias("dP0","tp1.fP[0]-tp2.fP[0]");
50   chainTr->SetAlias("dP1","tp1.fP[1]-tp2.fP[1]");
51   chainTr->SetAlias("dP2","tp1.fP[2]-tp2.fP[2]");
52   chainTr->SetAlias("dP3","tp1.fP[3]-tp2.fP[3]");
53   chainTr->SetAlias("dP4","tp1.fP[4]-tp2.fP[4]");
54   //
55   chainTr->SetAlias("sP0","sqrt(tp1.fC[0]+tp2.fC[0])");
56   chainTr->SetAlias("sP1","sqrt(tp1.fC[2]+tp2.fC[2])");
57   chainTr->SetAlias("sP2","sqrt(tp1.fC[5]+tp2.fC[5])");
58   chainTr->SetAlias("sP3","sqrt(tp1.fC[9]+tp2.fC[9])");
59   chainTr->SetAlias("sP4","sqrt(tp1.fC[14]+tp2.fC[14])");
60   //
61   chainTr->SetAlias("pP0","dP0/sP0");
62   chainTr->SetAlias("pP1","dP1/sP1");
63   chainTr->SetAlias("pP2","dP2/sP2");
64   chainTr->SetAlias("pP3","dP3/sP3");
65   //
66   chainTr->SetAlias("side","(sign(tp1.fP[1])+0)");
67   chainTr->SetAlias("dR","(1-abs(tp1.fP[1]/250.))");
68   chainTr->SetAlias("ta0","(tp1.fP[2]+tp2.fP[2])*0.5");
69   chainTr->SetAlias("ta1","(tp1.fP[3]+tp2.fP[3])*0.5");
70   chainTr->SetAlias("ca","cos(tp1.fAlpha+0)");
71   chainTr->SetAlias("sa","sin(tp1.fAlpha+0)");
72   //
73   chainTr->SetAlias("meanZ","(tp1.fP[1]+tp2.fP[1])*0.5");
74   chainTr->SetAlias("vx1","(v1.fElements[0]+0)");
75   chainTr->SetAlias("vy1","(v1.fElements[1]+0)");
76   chainTr->SetAlias("vz1","(v1.fElements[2]+0)");
77   chainTr->SetAlias("vdydx1","(v1.fElements[3]+0)");
78   chainTr->SetAlias("vdzdx1","(v1.fElements[4]+0)");
79   chainTr->SetAlias("vx2","(v2.fElements[0]+0)");
80   chainTr->SetAlias("vy2","(v2.fElements[1]+0)");
81   chainTr->SetAlias("vz2","(v2.fElements[2]+0)");
82   chainTr->SetAlias("vdydx2","(v2.fElements[3]+0)");
83   chainTr->SetAlias("vdzdx2","(v2.fElements[4]+0)");
84   //
85   //
86   //
87   chainTr->SetAlias("dx1","(AliTPCcalibAlign::SCorrect(0,0,s1,s2,vx1,vy1,vz1,vdydx1,vdzdx1)-vx1)");
88   chainTr->SetAlias("dy1","(AliTPCcalibAlign::SCorrect(0,1,s1,s2,vx1,vy1,vz1,vdydx1,vdzdx1)-vy1)");
89   chainTr->SetAlias("dz1","(AliTPCcalibAlign::SCorrect(0,2,s1,s2,vx1,vy1,vz1,vdydx1,vdzdx1)-vz1)");
90   chainTr->SetAlias("ddy1","(AliTPCcalibAlign::SCorrect(0,3,s1,s2,vx1,vy1,vz1,vdydx1,vdzdx1)-vdydx1)");
91   chainTr->SetAlias("ddz1","(AliTPCcalibAlign::SCorrect(0,4,s1,s2,vx1,vy1,vz1,vdydx1,vdzdx1)-vdzdx1)");
92   chainTr->SetAlias("ddy01","(AliTPCcalibAlign::SCorrect(0,5,s1,s2,vx1,vy1,vz1,vdydx1,vdzdx1)-vdydx1)");
93   chainTr->SetAlias("ddz01","(AliTPCcalibAlign::SCorrect(0,6,s1,s2,vx1,vy1,vz1,vdydx1,vdzdx1)-vdzdx1)");
94
95
96 }
97
98
99 void InitCutsAlign(){
100   /// resolution cuts
101
102   TCut cutS0("sqrt(tp2.fC[0]+tp2.fC[0])<0.2");
103   TCut cutS1("sqrt(tp2.fC[2]+tp2.fC[2])<0.2");
104   TCut cutS2("sqrt(tp2.fC[5]+tp2.fC[5])<0.01");
105   TCut cutS3("sqrt(tp2.fC[9]+tp2.fC[9])<0.01");
106   TCut cutS4("sqrt(tp2.fC[14]+tp2.fC[14])<0.5");
107   TCut cutS=cutS0+cutS1+cutS2+cutS3+cutS4;
108   //
109   // parameters matching cuts
110   TCut cutP0("abs(tp1.fP[0]-tp2.fP[0])<0.6");
111   TCut cutP1("abs(tp1.fP[1]-tp2.fP[1])<0.6");
112   TCut cutP2("abs(tp1.fP[2]-tp2.fP[2])<0.03");
113   TCut cutP3("abs(tp1.fP[3]-tp2.fP[3])<0.03");
114   TCut cutP=cutP0+cutP1+cutP2+cutP3;
115   //
116   TCut cutAll = cutS+cutP;
117   chainTr->Draw(">>listELP",cutAll,"entryList");
118   TEntryList *elist = (TEntryList*)gDirectory->Get("listELP");
119   chainTr->SetEntryList(elist);
120   
121   TCut cutRun("1");
122   TCut cutN120("1");
123
124 }
125
126 void MakeAlign(){
127   
128   align.ProcessTree(chainTr);
129   align.EvalFitters();
130   align.MakeTree("alignTree.root");
131   align.SetInstance(&align);
132 }
133  
134 void MakeCompareAlign(){
135   ///
136
137   TFile falignTreeNoMag("/lustre_alpha/alice/miranov/rec/LHC08d/nomag/alignTree.root");
138   TTree * treeAlignNoMag = (TTree*)falignTreeNoMag.Get("Align");
139   TFile falignTree("alignTree.root");
140   TTree * treeAlign = (TTree*)falignTree.Get("Align");
141   treeAlign->AddFriend(treeAlignNoMag,"NoMag");
142   treeAlignNoMag->SetMarkerStyle(26);
143   treeAlign->SetMarkerStyle(25);
144 }
145
146 TMatrixD * arrayAlign[72];
147 TMatrixD * arrayAlignTmp[72];
148
149 void ClearMatrices(){
150   ///
151
152   for (Int_t i=0;i<72; i++) {
153     TMatrixD * mat = new TMatrixD(4,4);
154     mat->UnitMatrix();
155     arrayAlign[i]=mat;
156     arrayAlignTmp[i]=(TMatrixD*)(mat->Clone());
157   }
158
159 }
160
161 void GlobalAlign(){
162   /// Global Align
163
164   TTreeSRedirector *cstream = new TTreeSRedirector("galign.root");
165
166   for (Int_t iter=0; iter<10;iter++){
167     printf("Iter=\t%d\n",iter);
168     for (Int_t is0=0;is0<72; is0++) {
169       //
170       TMatrixD  *mati0 = arrayAlign[is0];
171       TMatrixD matDiff(4,4);
172       Double_t sumw=0;      
173       for (Int_t is1=0;is1<72; is1++) {
174         Bool_t invers=kFALSE;
175         const TMatrixD *mat = align.GetTransformation(is0,is1,0); 
176         if (!mat){
177           invers=kTRUE;
178           mat = align.GetTransformation(is1,is0,0); 
179         }
180         if (!mat) continue;
181         Double_t weight=1;
182         if  ( (is1%18-is0%18)!=0) weight*=0.3;
183         if (is1/36>is0/36) weight*=2./3.; 
184         if (is1/36<is0/36) weight*=1./3.;
185         //
186         //
187         TMatrixD matT = *mat;   
188         if (invers)  matT.Invert();
189         matDiff+=weight*(*(arrayAlign[is1]))*matT;
190         sumw+=weight;
191         (*cstream)<<"LAlign"<<
192           "iter="<<iter<<
193           "s0="<<is0<<
194           "m6.="<<arrayAlign[is0]<<
195           "\n";
196
197       }
198       if (sumw>0){
199         matDiff*=1/sumw;
200         (*arrayAlignTmp[is0]) = matDiff;
201       }
202     }
203     for (Int_t is0=0;is0<72; is0++) {
204       (*arrayAlign[is0]) = (*arrayAlignTmp[is0]);
205       //      TMatrixD * matM1=  align.GetTransformation(is0,36+(is0+35)%36,0);
206       //TMatrixD * mat  =  align.GetTransformation(is0,36+(is0+36)%36,0);
207       //TMatrixD * matP1=  align.GetTransformation(is0,36+(is0+37)%36,0);
208       //
209       (*cstream)<<"GAlign"<<
210         "iter="<<iter<<
211         "s0="<<is0<<
212         "m6.="<<arrayAlign[is0]<<
213         "\n";
214     }    
215   }
216   delete cstream;
217 }
218
219
220
221
222
223 void MakeGlobalCorr(){
224   ///
225
226   TStatToolkit toolkit;
227   Double_t chi2=0;
228   Int_t    npoints=0;
229   TVectorD fitParam;
230   TMatrixD covMatrix;
231   TVectorD chi2V(5);
232   //
233   TString fstring="";
234   fstring+="side++";
235   //
236   fstring+="dR++";
237   fstring+="dR*dR++";
238   fstring+="dR*sa++";
239   fstring+="dR*ca++";
240   fstring+="ta0++";
241   fstring+="ta1++";
242   //
243   fstring+="dR*side++";
244   fstring+="dR*dR*side++";
245   fstring+="dR*sa*side++";
246   fstring+="dR*ca*side++";
247   fstring+="ta0*side++";
248   fstring+="ta1*side++";
249
250   TString * strP0 = TStatToolkit::FitPlane(chainTr,"dP0:sP0", fstring.Data(),cutS, chi2,npoints,fitParam,covMatrix);
251   chi2V[0]=TMath::Sqrt(chi2/npoints);
252   chainTr->SetAlias("corrP0",strP0->Data());
253
254   TString * strP1 = TStatToolkit::FitPlane(chainTr,"dP1:sP1", fstring.Data(),cutS, chi2,npoints,fitParam,covMatrix);
255   chi2V[1]=TMath::Sqrt(chi2/npoints);
256   chainTr->SetAlias("corrP1",strP1->Data());
257
258   TString * strP2 = TStatToolkit::FitPlane(chainTr,"dP2:sP2", fstring.Data(),cutS, chi2,npoints,fitParam,covMatrix);
259   chi2V[2]=TMath::Sqrt(chi2/npoints);
260   chainTr->SetAlias("corrP2",strP2->Data());
261
262   TString * strP3 = TStatToolkit::FitPlane(chainTr,"dP3:sP3", fstring.Data(),cutS, chi2,npoints,fitParam,covMatrix);
263   chi2V[3]=TMath::Sqrt(chi2/npoints);
264   chainTr->SetAlias("corrP3",strP3->Data());
265
266   TString * strP4 = TStatToolkit::FitPlane(chainTr,"dP4:sP4", fstring.Data(),cutS, chi2,npoints,fitParam,covMatrix);
267   chi2V[4]=TMath::Sqrt(chi2/npoints);
268   chainTr->SetAlias("corrP4",strP4->Data());
269 }
270
271 void P0resolZ(){
272   ///
273
274   TH2F * hdP0Z = new TH2F("hdP0Z","hdP0Z",10,-250,250,100,-0.5,0.5);
275   TH2F * hdP0ZNoCor = new TH2F("hdP0ZNoCor","hdP0ZNoCor",10,-250,250,100,-0.5,0.5);
276   chainTr->Draw("(dP0-corrP0)/sqrt(2.):meanZ>>hdP0Z",""+cutRun+cutS+cutN120,"");
277   chainTr->Draw("(dP0-0)/sqrt(2.):meanZ>>hdP0ZNoCor",""+cutRun+cutS+cutN120,"");
278   
279   hdP0Z->FitSlicesY();
280   hdP0ZNoCor->FitSlicesY();
281   hdP0Z_2->SetMinimum(0);
282   hdP0Z_2->SetXTitle("Z position (cm)");
283   hdP0Z_2->SetYTitle("#sigma_{y} (cm)");
284   hdP0Z_2->SetMarkerStyle(25);
285   hdP0ZNoCor_2->SetMarkerStyle(26);
286   hdP0Z_2->Draw();
287   hdP0ZNoCor_2->Draw("same");  
288   gPad->SaveAs("picAlign/SigmaY_z.gif");
289   gPad->SaveAs("picAlign/SigmaY_z.eps");
290   //
291   hdP0ZNoCor_1->SetXTitle("Z position (cm)");
292   hdP0ZNoCor_1->SetYTitle("#Delta{y} (cm)");
293   hdP0Z_1->SetMarkerStyle(25);
294   hdP0ZNoCor_1->SetMarkerStyle(26);
295   hdP0ZNoCor_1->Draw("");  
296   hdP0Z_1->Draw("same");
297   gPad->SaveAs("picAlign/DeltaY_z.gif");
298   gPad->SaveAs("picAlign/DeltaY_z.eps");
299   //
300   //
301   TH2F * hdPP0Z = new TH2F("hdPP0Z","hdPP0Z",8,-200,200,50,-5.05,5.05);
302   TH2F * hncdPP0Z = new TH2F("hncdPP0Z","hncdPP0Z",8,-200,200,50,-5.05,5.05);
303   chainTr->Draw("(dP0-corrP0)/sP0:meanZ>>hdPP0Z",""+cutRun+cutS+cutN120,"");
304   chainTr->Draw("(dP0-0)/sP0:meanZ>>hncdPP0Z",""+cutRun+cutS+cutN120,"");
305   hdPP0Z->FitSlicesY();
306   hncdPP0Z->FitSlicesY();
307   hdPP0Z_2->SetMarkerStyle(25);
308   hncdPP0Z_2->SetMarkerStyle(26);
309   hdPP0Z_2->SetMinimum(0);
310   hdPP0Z_2->SetXTitle("Z position (cm)");
311   hdPP0Z_2->SetYTitle("#sigma y (Unit) ");
312   hdPP0Z_2->Draw();
313   hncdPP0Z_2->Draw("same");
314   gPad->SaveAs("picAlign/PullY_z.gif");
315   gPad->SaveAs("picAlign/PullY_z.eps");
316 }
317
318 void P1resolZ(){
319   ///
320
321   TH2F * hdP1Z = new TH2F("hdP1Z","hdP1Z",10,-250,250,100,-0.2,0.2);
322   TH2F * hdP1ZNoCor = new TH2F("hdP1ZNoCor","hdP1ZNoCor",10,-250,250,100,-0.2,0.2);
323   chainTr->Draw("(dP1-corrP1)/sqrt(2.):meanZ>>hdP1Z",""+cutRun+cutS+cutN120,"");
324   chainTr->Draw("(dP1-0)/sqrt(2.):meanZ>>hdP1ZNoCor",""+cutRun+cutS+cutN120,"");
325   
326   hdP1Z->FitSlicesY();
327   hdP1ZNoCor->FitSlicesY();
328   hdP1Z_2->SetMinimum(0);
329   hdP1Z_2->SetXTitle("Z position (cm)");
330   hdP1Z_2->SetYTitle("#sigma_{z} (cm)");
331   hdP1Z_2->SetMarkerStyle(25);
332   hdP1ZNoCor_2->SetMarkerStyle(26);
333   hdP1Z_2->Draw();
334   hdP1ZNoCor_2->Draw("same");  
335   gPad->SaveAs("picAlign/SigmaZ_z.gif");
336   gPad->SaveAs("picAlign/SigmaZ_z.eps");
337   //
338   hdP1ZNoCor_1->SetXTitle("Z position (cm)");
339   hdP1ZNoCor_1->SetYTitle("#Delta{z} (cm)");
340   hdP1Z_1->SetMarkerStyle(25);
341   hdP1ZNoCor_1->SetMarkerStyle(26);
342   hdP1ZNoCor_1->Draw("");  
343   hdP1Z_1->Draw("same");
344   gPad->SaveAs("picAlign/DeltaZ_z.gif");
345   gPad->SaveAs("picAlign/DeltaZ_z.eps");
346   //
347   //
348   TH2F * hdPP1Z = new TH2F("hdPP1Z","hdPP1Z",8,-200,200,50,-5.05,5.05);
349   TH2F * hncdPP1Z = new TH2F("hncdPP1Z","hncdPP1Z",8,-200,200,50,-5.05,5.05);
350   chainTr->Draw("(dP1-corrP1)/sP1:meanZ>>hdPP1Z",""+cutRun+cutS+cutN120,"");
351   chainTr->Draw("(dP1-0)/sP1:meanZ>>hncdPP1Z",""+cutRun+cutS+cutN120,"");
352   hdPP1Z->FitSlicesY();
353   hncdPP1Z->FitSlicesY();
354   hdPP1Z_2->SetMarkerStyle(25);
355   hncdPP1Z_2->SetMarkerStyle(26);
356   hdPP1Z_2->SetMinimum(0);
357   hdPP1Z_2->SetXTitle("Z position (cm)");
358   hdPP1Z_2->SetYTitle("#sigma z (Unit) ");
359   hdPP1Z_2->Draw();
360   hncdPP1Z_2->Draw("same");
361   gPad->SaveAs("picAlign/PullZ_z.gif");
362   gPad->SaveAs("picAlign/PullZ_z.eps");
363 }
364
365
366 void P4resolZ(){
367   ///
368
369   TH2F * hdP4Z = new TH2F("hdP4Z","hdP4Z",10,-250,250,100,-0.4,0.4);
370   TH2F * hdP4ZNoCor = new TH2F("hdP4ZNoCor","hdP4ZNoCor",10,-250,250,100,-0.4,0.4);
371   chainTr->Draw("(dP4-corrP4)/sqrt(2.):meanZ>>hdP4Z",""+cutRun+cutS+cutN120,"");
372   chainTr->Draw("(dP4-0)/sqrt(2.):meanZ>>hdP4ZNoCor",""+cutRun+cutS+cutN120,"");
373   
374   hdP4Z->FitSlicesY();
375   hdP4ZNoCor->FitSlicesY();
376   hdP4Z_2->SetMinimum(0);
377   hdP4Z_2->SetXTitle("Z position (cm)");
378   hdP4Z_2->SetYTitle("#sigma_{1/pt} (1/GeV)");
379   hdP4Z_2->SetMarkerStyle(25);
380   hdP4ZNoCor_2->SetMarkerStyle(26);
381   hdP4Z_2->Draw();
382   hdP4ZNoCor_2->Draw("same");  
383   gPad->SaveAs("picAlign/SigmaP4_z.gif");
384   gPad->SaveAs("picAlign/SigmaP4_z.eps");
385   //
386   hdP4ZNoCor_1->SetXTitle("Z position (cm)");
387   hdP4ZNoCor_1->SetYTitle("#Delta_{1/p_{t}} (1/GeV)");
388   hdP4Z_1->SetMarkerStyle(25);
389   hdP4ZNoCor_1->SetMarkerStyle(26);
390   hdP4ZNoCor_1->Draw("");  
391   hdP4Z_1->Draw("same");
392   gPad->SaveAs("picAlign/DeltaP4_z.gif");
393   gPad->SaveAs("picAlign/DeltaP4_z.eps");
394   //
395   //
396   TH2F * hdPP4Z = new TH2F("hdPP4Z","hdPP4Z",8,-200,200,50,-5.05,5.05);
397   TH2F * hncdPP4Z = new TH2F("hncdPP4Z","hncdPP4Z",8,-200,200,50,-5.05,5.05);
398   chainTr->Draw("(dP4-corrP4)/sP4:meanZ>>hdPP4Z",""+cutRun+cutS+cutN120,"");
399   chainTr->Draw("(dP4-0)/sP4:meanZ>>hncdPP4Z",""+cutRun+cutS+cutN120,"");
400   hdPP4Z->FitSlicesY();
401   hncdPP4Z->FitSlicesY();
402   hdPP4Z_2->SetMarkerStyle(25);
403   hncdPP4Z_2->SetMarkerStyle(26);
404   hdPP4Z_2->SetMinimum(0);
405   hdPP4Z_2->SetXTitle("Z position (cm)");
406   hdPP4Z_2->SetYTitle("#sigma 1/p_{t} (Unit) ");
407   hdPP4Z_2->Draw();
408   hncdPP4Z_2->Draw("same");
409   gPad->SaveAs("picAlign/PullP4_z.gif");
410   gPad->SaveAs("picAlign/PullP4_z.eps");
411 }
412
413
414
415
416
417 void MakeFit(){
418   ///
419
420   TChain *chainTracklet=chainTr;
421   AliTPCcalibAlign align;
422   //
423   TVectorD * vec1 = 0;
424   TVectorD * vec2 = 0;
425   AliExternalTrackParam * tp1 = 0;
426   AliExternalTrackParam * tp2 = 0;  
427   Int_t      s1 = 0;
428   Int_t      s2 = 0;
429   {
430     for (Int_t i=0; i< elist->GetN(); i++){
431       //for (Int_t i=0; i<100000; i++){
432       chainTracklet->GetBranch("tp1.")->SetAddress(&tp1);
433       chainTracklet->GetBranch("tp2.")->SetAddress(&tp2);
434       chainTracklet->GetBranch("v1.")->SetAddress(&vec1);
435       chainTracklet->GetBranch("v2.")->SetAddress(&vec2);
436       chainTracklet->GetBranch("s1")->SetAddress(&s1);
437       chainTracklet->GetBranch("s2")->SetAddress(&s2);      
438       chainTracklet->GetEntry(i);
439       if (i%100==0) printf("%d\t%d\t%d\t\n",i, s1,s2);
440       //vec1.Print();
441       TLinearFitter * fitter = align.GetOrMakeFitter6(s1,s2);
442       if (fitter) align.Process6(vec1->GetMatrixArray(),vec2->GetMatrixArray(),fitter);
443     }
444   }
445 }
446
447
448 // chainTr->Scan("vy1:AliTPCcalibAlign::SCorrect(0,1,s1,s2,vx1,vy1,vz1,vdydx1,vdzdx1)","","",100); 
449
450 void MakePlot(){
451   ///
452
453   chainTr->Draw("vy2-vy1:s1>>hisdYNoCor(36,0,36,100,-0.3,0.3)","abs(s2-s1-36)<1","",40000);
454   chainTr->Draw("vy2-vy1-dy1:s1>>hisdYCor(36,0,36,100,-0.3,0.3)","abs(s2-s1-36)<1","",40000);
455   hisdYCor->FitSlicesY();
456   hisdYNoCor->FitSlicesY();
457   hisdYCor_1->Draw("");
458
459 }
460
461
462 void dPhi(){
463   ///
464
465   treeAlign->Draw("1000*m6.fElements[4]","s2==s1+36&&nphi>100");
466   htemp->SetXTitle("#Delta_{#phi} (mrad)");
467   gPad->SaveAs("picAlign/mag5dPhi.eps");
468   gPad->SaveAs("picAlign/mag5dPhi.gif");
469   //
470   treeAlignNoMag->Draw("1000*m6.fElements[4]","s2==s1+36&&nphi>100");
471   htemp->SetXTitle("#Delta_{#phi} (mrad)");
472   gPad->SaveAs("picAlign/nomagdPhi.eps");
473   gPad->SaveAs("picAlign/nomagdPhi.gif");
474   //
475   //
476   //
477   treeAlign->Draw("1000*(m6.fElements[4]-NoMag.m6.fElements[4])","s2==s1+36&&nphi>100");
478   htemp->SetXTitle("#Delta_{#phi} (mrad)");
479   gPad->SaveAs("picAlign/diffnomagmag5dPhi.eps");
480   gPad->SaveAs("picAlign/diffnomagmag5dPhi.gif");
481   //
482 }
483
484
485 void dTheta(){
486   ///
487
488   treeAlign->Draw("1000*m6.fElements[8]","s2==s1+36&&nphi>100");
489   htemp->SetXTitle("#Delta_{#theta} (mrad)");
490   gPad->SaveAs("picAlign/mag5dTheta.eps");
491   gPad->SaveAs("picAlign/mag5dTheta.gif");
492   //
493   treeAlignNoMag->Draw("1000*m6.fElements[8]","s2==s1+36&&nphi>100");
494   htemp->SetXTitle("#Delta_{#theta} (mrad)");
495   gPad->SaveAs("picAlign/nomagdTheta.eps");
496   gPad->SaveAs("picAlign/nomagdTheta.gif");
497   //
498   //
499   //
500   treeAlign->Draw("1000*(m6.fElements[8]-NoMag.m6.fElements[8])","s2==s1+36&&nphi>100");
501   htemp->SetXTitle("#Delta_{#theta} (mrad)");
502   gPad->SaveAs("picAlign/diffnomagmag5dTheta.eps");
503   gPad->SaveAs("picAlign/diffnomagmag5dTheta.gif");
504   //
505   //
506   treeAlign->Draw("1000*(m6.fElements[8]):s1","s2==s1+36&&nphi>100");
507   htemp->SetYTitle("#Delta_{#theta} (mrad)");
508   htemp->SetXTitle("Sector number");
509   treeAlignNoMag->Draw("1000*(m6.fElements[8]):s1","s2==s1+36&&nphi>100","same");
510   gPad->SaveAs("picAlign/diffnomagmag5dTheta_s1.eps");
511   gPad->SaveAs("picAlign/diffnomagmag5dTheta_s1.gif");
512   //
513 }
514
515
516
517 void dZ(){
518   ///
519
520   treeAlign->Draw("dz","s2==s1+36&&nphi>100");
521   htemp->SetXTitle("#Delta_{Z} (cm)");
522   gPad->SaveAs("picAlign/mag5dZ.eps");
523   gPad->SaveAs("picAlign/mag5dZ.gif");
524   //
525   treeAlignNoMag->Draw("dz","s2==s1+36&&nphi>100");
526   htemp->SetXTitle("#Delta_{Z} (cm)");
527   gPad->SaveAs("picAlign/nomagdZ.eps");
528   gPad->SaveAs("picAlign/nomagdZ.gif");
529   //
530   //
531   //
532   treeAlign->Draw("(dz-NoMag.dz)","s2==s1+36&&nphi>100");
533   htemp->SetXTitle("#Delta_{Z} (cm)");
534   gPad->SaveAs("picAlign/diffnomagmag5dZ.eps");
535   gPad->SaveAs("picAlign/diffnomagmag5dZ.gif");
536   //
537   //
538   treeAlign->Draw("dz:s1","s2==s1+36&&nphi>100");
539   htemp->SetYTitle("#Delta_{Z} (cm)");
540   htemp->SetXTitle("Sector number");
541   treeAlignNoMag->Draw("dz:s1","s2==s1+36&&nphi>100","same");
542   gPad->SaveAs("picAlign/diffnomagmag5dZ_s1.eps");
543   gPad->SaveAs("picAlign/diffnomagmag5dZ_s1.gif");
544   //
545 }
546
547 void dY(){
548   ///
549
550   treeAlign->Draw("dy","s2==s1+36&&nphi>100");
551   htemp->SetXTitle("#Delta_{Y} (cm)");
552   gPad->SaveAs("picAlign/mag5dY.eps");
553   gPad->SaveAs("picAlign/mag5dY.gif");
554   //
555   treeAlignNoMag->Draw("dy","s2==s1+36&&nphi>100");
556   htemp->SetXTitle("#Delta_{Y} (cm)");
557   gPad->SaveAs("picAlign/nomagdY.eps");
558   gPad->SaveAs("picAlign/nomagdY.gif");
559   //
560   //
561   //
562   treeAlign->Draw("(dy-NoMag.dy)","s2==s1+36&&nphi>100");
563   htemp->SetXTitle("#Delta_{Y} (cm)");
564   gPad->SaveAs("picAlign/diffnomagmag5dY.eps");
565   gPad->SaveAs("picAlign/diffnomagmag5dY.gif");
566   //
567   //
568   treeAlign->Draw("dy:s1","s2==s1+36&&nphi>100");
569   htemp->SetYTitle("#Delta_{Y} (cm)");
570   htemp->SetXTitle("Sector number");
571   treeAlignNoMag->Draw("dy:s1","s2==s1+36&&nphi>100","same");
572   gPad->SaveAs("picAlign/diffnomagmag5dY_s1.eps");
573   gPad->SaveAs("picAlign/diffnomagmag5dY_s1.gif");
574   //
575 }
576
577
578
579
580
581
582
583
584