]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/macrosSDD/CheckCalibs.C
Overlaps corrected, new shape of sectors
[u/mrichter/AliRoot.git] / ITS / macrosSDD / CheckCalibs.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TApplication.h>
3 #include <TGClient.h>
4 #include <TGButton.h>
5 #include <TGFrame.h>
6 #include <TGLayout.h>
7 #include <TGWindow.h>
8 #include <TGLabel.h>
9 #include <TGNumberEntry.h>
10 #include <TString.h>
11 #include <TCanvas.h>
12 #include <TStyle.h>
13 #include <TPad.h>
14 #include <TH1F.h>
15 #include <TH2F.h>
16 #include <TF1.h>
17 #include <TSystem.h>
18 #include <TGraph.h>
19 #include <TLatex.h>
20 #include <TLine.h>
21 #include <TROOT.h>
22 #endif
23
24 // Macro to display the output of SDD calibration runs
25 // -> copies the ASCII files with the DA output from the LDC
26 // -> plot various histograms and graphs
27 // Origin: F. Prino (prino@to.infn.it)
28
29 class CheckCalibInterface : public TGMainFrame {
30
31 private:
32   TGCompositeFrame    *fHor1;
33   TGCompositeFrame    *fHor2;
34   TGCompositeFrame    *fHor3;
35   TGCompositeFrame    *fHor4;
36   TGTextButton        *fCopyFiles;
37   TGTextButton        *fShowPedest;
38   TGTextButton        *fShowPulser;
39   TGTextButton        *fShowInject;
40   TGTextButton        *fShowOneMod;
41   TGTextButton        *fShowDDL;
42   TGTextButton        *fExit;
43   TGGroupFrame        *fGframeALL;
44   TGGroupFrame        *fGframeSING;
45   TGGroupFrame        *fGframe1;
46   TGGroupFrame        *fGframe2;
47   TGNumberEntry       *fDDL;
48   TGNumberEntry       *fChannel;
49   Int_t fNumDDL;
50   Int_t fNumChannel;
51
52 public:
53   CheckCalibInterface(const TGWindow *p, UInt_t w, UInt_t h);
54   virtual ~CheckCalibInterface();
55   void DoSetlabel();
56
57   static void CopyFiles();
58   static void ShowPedestal();
59   static void ShowPulser();
60   static void ShowInjector();
61   static void ShowSingleModule(Int_t iddl, Int_t ichan);
62   static void ShowDDL(Int_t iddl);
63   static void ClearAll();
64
65   ClassDef(CheckCalibInterface, 0)
66 };
67                           
68 CheckCalibInterface::CheckCalibInterface(const TGWindow *p, UInt_t w, UInt_t h)
69    : TGMainFrame(p, w, h)
70 {
71
72   fHor1 = new TGHorizontalFrame(this);
73   fHor2 = new TGHorizontalFrame(this);
74   fHor3 = new TGHorizontalFrame(this);
75   fHor4 = new TGHorizontalFrame(this);
76
77   TGLayoutHints *lh=new TGLayoutHints(kLHintsTop | kLHintsLeft, 5, 5, 5, 5);
78   TGLayoutHints *lhc=new TGLayoutHints(kLHintsCenterY | kLHintsExpandY, 5, 5, 5, 5);
79   TGLayoutHints *lh1=new TGLayoutHints(kLHintsTop | kLHintsExpandX, 5, 5, 5, 5);
80
81   fCopyFiles = new TGTextButton(fHor1, "&Copy Files", "CheckCalibInterface::CopyFiles()");
82   fHor1->AddFrame(fCopyFiles, lh1);
83
84   fGframeALL = new TGGroupFrame(fHor2,"All Modules",kHorizontalFrame);
85   fShowPedest = new TGTextButton(fGframeALL, "&Show Pedestal","CheckCalibInterface::ShowPedestal()");
86   fShowPulser = new TGTextButton(fGframeALL, "&Show Pulser","CheckCalibInterface::ShowPulser()");
87   fShowInject = new TGTextButton(fGframeALL, "&Show Injector","CheckCalibInterface::ShowInjector()");    
88   fGframeALL->AddFrame(fShowPedest, lh1);
89   fGframeALL->AddFrame(fShowPulser, lh1);
90   fGframeALL->AddFrame(fShowInject, lh1);
91   fHor2->AddFrame(fGframeALL, lh1);
92
93   fGframeSING = new TGGroupFrame(fHor3,"Single Module",kHorizontalFrame);
94   fGframe1 = new TGGroupFrame(fGframeSING, "DDL");
95   fDDL = new TGNumberEntry(fGframe1, 0, 2,-1, TGNumberFormat::kNESInteger,
96                            TGNumberFormat::kNEANonNegative, 
97                            TGNumberFormat::kNELLimitMinMax,
98                            0, 23);
99   fGframe2 = new TGGroupFrame(fGframeSING, "Channel");
100   fChannel = new TGNumberEntry(fGframe2, 0, 2,-1, TGNumberFormat::kNESInteger,
101                                TGNumberFormat::kNEANonNegative, 
102                                TGNumberFormat::kNELLimitMinMax,
103                                0, 11);
104   fGframe1->AddFrame(fDDL,lh);
105   fGframe2->AddFrame(fChannel,lh);
106   fGframeSING->AddFrame(fGframe1,lh);
107   fGframeSING->AddFrame(fGframe2,lh);
108   
109   fDDL->Connect("ValueSet(Long_t)", "CheckCalibInterface", this, "DoSetlabel()");
110   (fDDL->GetNumberEntry())->Connect("ReturnPressed()", "CheckCalibInterface", this, "DoSetlabel()");
111   
112   fChannel->Connect("ValueSet(Long_t)", "CheckCalibInterface", this, "DoSetlabel()");
113   (fChannel->GetNumberEntry())->Connect("ReturnPressed()", "CheckCalibInterface", this, "DoSetlabel()");
114   fHor3->AddFrame(fGframeSING, lh1);
115   
116   fShowOneMod = new TGTextButton(fGframeSING, "&Show Selected Module","CheckCalibInterface::ShowSingleModule(0,0)");
117   fGframeSING->AddFrame(fShowOneMod, lhc);
118
119   fShowDDL = new TGTextButton(fGframeSING, "&Show vdrift for DDL","CheckCalibInterface::ShowDDL(0)");
120   fGframeSING->AddFrame(fShowDDL, lhc);
121
122   fExit = new TGTextButton(fHor4, "&Exit", "gApplication->Terminate(0)");
123   fHor4->AddFrame(fExit, lh1);
124
125   AddFrame(fHor1,lh1);
126   AddFrame(fHor2,lh1);
127   AddFrame(fHor3,lh1);
128   AddFrame(fHor4,lh1);
129
130   SetCleanup(kDeepCleanup);
131   SetWindowName("Main Control");
132   MapSubwindows();
133   Resize(GetDefaultSize());
134   MapWindow();
135 }
136
137 CheckCalibInterface::~CheckCalibInterface()
138 {
139    // Destructor.
140    
141    Cleanup();
142 }
143
144
145 void CheckCalibInterface::CopyFiles(){
146   //
147   // 8 LDCs in Run2 with new (simpler) dirs - M.S. 25 jul 14
148   TString command;
149   TString ldcName[8]={"aldaqpc009","aldaqpc010","aldaqpc011","aldaqpc012",
150                       "aldaqpc013","aldaqpc014","aldaqpc015","aldaqpc016"};
151   gSystem->Exec("rm -rf calibFiles");
152   gSystem->Exec("mkdir calibFiles");
153   for(Int_t iLDC=0; iLDC<8; iLDC++){
154     Int_t firstDDL=iLDC*3;
155     Int_t lastDDL=iLDC*3+2;
156     command.Form("scp %s:/dateSite/ldc-SDD-%d/work/SDDbase_step2_LDC.tar calibFiles/.",ldcName[iLDC].Data(),iLDC);
157     gSystem->Exec(command.Data());
158     command.Form("scp %s:/dateSite/ldc-SDD-%d/work/SDDbase_LDC.tar calibFiles/.",ldcName[iLDC].Data(),iLDC);
159     gSystem->Exec(command.Data());
160     command.Form("scp %s:/dateSite/ldc-SDD-%d/work/SDDinj_LDC.tar calibFiles/.",ldcName[iLDC].Data(),iLDC);
161     gSystem->Exec(command.Data());
162     gSystem->Exec("cd calibFiles; tar xvf SDDbase_step2_LDC.tar; tar xvf SDDbase_LDC.tar; tar xvf SDDinj_LDC.tar; cd ..");
163   }
164   printf("-------------- DONE ---------------\n");
165
166   return;
167 }
168
169 void CheckCalibInterface::ShowPedestal(){
170   ClearAll();
171
172   TH2F* hbadchan=new TH2F("hbadchan","Number of bad channels",24,-0.25,11.75,24,-0.5,23.5);
173   TH2F* hrawnoisemod=new TH2F("hrawnoisemod","",288,-0.5,287.5,50,0.,10.);
174   TH1F* hrawnoise=new TH1F("hrawnoise","",100,0.,10.);
175   TH2F* hcornoisemod=new TH2F("hcornoisemod","",288,-0.5,287.5,50,0.,10.);
176   TH1F* hcornoise=new TH1F("hcornoise","",100,0.,10.);
177   TH2F* hbasemod=new TH2F("hbasemod","",288,-0.5,287.5,50,0.,150.);
178   TH1F* hbase=new TH1F("hbase","",100,0.,150.);
179   Int_t retfscf;
180   TString inpFileName;
181   Float_t baseline,rawnoise,cmn,corn;
182   Int_t isgoodan,i,basmin,basoff;
183   Int_t th,tl;
184   Int_t nGoodAnodes=0;
185   for(Int_t iddl=0;iddl<24;iddl++){
186     for(Int_t imod=0;imod<12;imod++){
187       for(Int_t isid=0;isid<2;isid++){
188         inpFileName.Form("./calibFiles/SDDbase_step2_ddl%02dc%02d_sid%d.data",iddl,imod,isid);
189         FILE* basFil = fopen(inpFileName.Data(),"read");
190         Int_t sideId=imod*2+isid;
191         Int_t modId=iddl*12+imod;
192         if (basFil == 0) hbadchan->SetBinContent(sideId+1,iddl+1,256);
193         else{
194           retfscf=fscanf(basFil,"%d\n",&th);
195           retfscf=fscanf(basFil,"%d\n",&tl);
196           if(th==255 && tl==20) continue;
197           for(Int_t ian=0;ian<256;ian++){
198             retfscf=fscanf(basFil,"%d %d %f %d %d %f %f %f\n",&i,&isgoodan,&baseline,&basmin,&basoff,&rawnoise,&cmn,&corn);
199             hrawnoisemod->Fill(modId,rawnoise);
200             hrawnoise->Fill(rawnoise);
201             hcornoisemod->Fill(modId,corn);
202             hcornoise->Fill(corn);
203             hbasemod->Fill(modId,baseline);
204             hbase->Fill(baseline);
205             if(!isgoodan){
206               hbadchan->SetBinContent(sideId+1,iddl+1, 1+hbadchan->GetBinContent(sideId+1,iddl+1));
207             }else{
208               nGoodAnodes++;
209             }
210           }
211           fclose(basFil);
212         }
213       }
214     }
215   }
216   hrawnoisemod->SetStats(0);
217   hcornoisemod->SetStats(0);
218   hbasemod->SetStats(0);
219   hbadchan->SetStats(0);
220   gStyle->SetPalette(1);
221   TString txtCountGood;
222   Float_t fracGood=100.*(Float_t)nGoodAnodes/(520.*260.);
223   txtCountGood.Form("Number of GoodAnodes = %d (%5.1f\%)\n",nGoodAnodes,fracGood);
224
225   TCanvas* c0=new TCanvas("c0","Bad Channels",800,900);
226   c0->SetRightMargin(0.14);
227   c0->SetBottomMargin(0.2);
228   hbadchan->Draw("colz"); 
229   hbadchan->GetXaxis()->SetTitle("Channel");
230   hbadchan->GetYaxis()->SetTitle("DDL");
231   hbadchan->GetXaxis()->SetTickLength(0);
232   hbadchan->GetYaxis()->SetTickLength(0);
233   c0->cd();
234   TLine** linv=new TLine*[12];
235   for(Int_t i=0;i<12;i++){
236     linv[i]=new TLine(i+0.75,-0.5,i+0.75,23.5);
237     linv[i]->SetLineColor(kGray+1);
238     linv[i]->Draw();
239   }
240   TLine** linh=new TLine*[24];
241   for(Int_t i=0;i<24;i++){
242     linh[i]=new TLine(-0.25,i+0.5,11.75,i+0.5);
243     linh[i]->SetLineColor(kGray+1);
244     linh[i]->Draw();
245   }
246   TLatex* tg=new TLatex(0.1,0.05,txtCountGood.Data());
247   tg->SetNDC();
248   tg->SetTextColor(4);
249   tg->SetTextSize(0.04);
250   tg->Draw();
251   c0->Update();
252  
253   TCanvas* c1=new TCanvas("c1","Baseline",1200,700);
254   c1->Divide(2,1);
255   c1->cd(1);
256   gPad->SetLeftMargin(0.14);
257   gPad->SetRightMargin(0.14);
258   hbase->Draw();
259   hbase->GetXaxis()->SetTitle("Baseline (ADC)");
260   c1->cd(2);
261   gPad->SetLeftMargin(0.14);
262   gPad->SetRightMargin(0.14);
263   hbasemod->Draw("colz");
264   hbasemod->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
265   hbasemod->GetYaxis()->SetTitle("Baseline (ADC)");
266   hbasemod->GetYaxis()->SetTitleOffset(1.35);
267
268   TCanvas* c2=new TCanvas("c2","Raw Noise",1200,700);
269   c2->Divide(2,1);
270   c2->cd(1);
271   gPad->SetLeftMargin(0.14);
272   gPad->SetRightMargin(0.14);
273   hrawnoise->Draw();
274   hrawnoise->GetXaxis()->SetTitle("Raw Noise (ADC)");
275   c2->cd(2);
276   gPad->SetLeftMargin(0.14);
277   gPad->SetRightMargin(0.14);
278   hrawnoisemod->Draw("colz");
279   hrawnoisemod->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
280   hrawnoisemod->GetYaxis()->SetTitle("Raw Noise (ADC)");
281   hrawnoisemod->GetYaxis()->SetTitleOffset(1.35);
282
283   TCanvas* c3=new TCanvas("c3","Corr Noise",1200,700);
284   c3->Divide(2,1);
285   c3->cd(1);
286   gPad->SetLeftMargin(0.14);
287   gPad->SetRightMargin(0.14);
288   hcornoise->Draw();
289   hcornoise->GetXaxis()->SetTitle("Raw Noise (ADC)");
290   c3->cd(2);
291   gPad->SetLeftMargin(0.14);
292   gPad->SetRightMargin(0.14);
293   hcornoisemod->Draw("colz");
294   hcornoisemod->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
295   hcornoisemod->GetYaxis()->SetTitle("Raw Noise (ADC)");
296   hcornoisemod->GetYaxis()->SetTitleOffset(1.35);  
297 }
298
299 void CheckCalibInterface::ShowPulser(){
300   ClearAll();
301
302   TH2F* hbadchan=new TH2F("hbadchan","Number of bad channels",24,-0.25,11.75,24,-0.5,23.5);
303   TH2F* hrawnoisemod=new TH2F("hrawnoisemod","",288,-0.5,287.5,50,0.,10.);
304   TH1F* hrawnoise=new TH1F("hrawnoise","",100,0.,10.);
305   TH2F* hcornoisemod=new TH2F("hcornoisemod","",288,-0.5,287.5,50,0.,10.);
306   TH1F* hcornoise=new TH1F("hcornoise","",100,0.,10.);
307   TH2F* hbasemod=new TH2F("hbasemod","",288,-0.5,287.5,50,0.,150.);
308   TH1F* hbase=new TH1F("hbase","",100,0.,150.);
309   TH2F* hgainmod=new TH2F("hgainmod","",288,-0.5,287.5,50,0.,5.);
310   TH1F* hgain=new TH1F("hgain","",100,0.,5.);
311
312   Int_t retfscf;
313   TString inpFileName;
314   Float_t baseline,rawnoise,cmn,corn,gain;
315   Int_t isgoodan,i,basmin,basoff;
316   Int_t th,tl;
317   Int_t iii,jjj,kkk;
318   Int_t nGoodAnodes=0;
319
320   for(Int_t iddl=0;iddl<24;iddl++){
321     for(Int_t imod=0;imod<12;imod++){
322       for(Int_t isid=0;isid<2;isid++){
323         inpFileName.Form("./calibFiles/SDDbase_ddl%02dc%02d_sid%d.data",iddl,imod,isid);
324         FILE* basFil = fopen(inpFileName.Data(),"read");
325         Int_t sideId=imod*2+isid;
326         Int_t modId=iddl*12+imod;
327         if (basFil == 0) hbadchan->SetBinContent(sideId+1,iddl+1,256);
328         else{
329           retfscf=fscanf(basFil,"%d %d %d\n",&iii,&jjj,&kkk);
330           if(kkk==0) continue;
331           retfscf=fscanf(basFil,"%d\n",&th);
332           retfscf=fscanf(basFil,"%d\n",&tl);
333           if(th==255 && tl==20) continue;
334           for(Int_t ian=0;ian<256;ian++){
335             retfscf=fscanf(basFil,"%d %d %f %d %d %f %f %f %f\n",&i,&isgoodan,&baseline,&basmin,&basoff,&rawnoise,&cmn,&corn,&gain);
336             hrawnoisemod->Fill(modId,rawnoise);
337             hrawnoise->Fill(rawnoise);
338             hcornoisemod->Fill(modId,corn);
339             hcornoise->Fill(corn);
340             hbasemod->Fill(modId,baseline);
341             hbase->Fill(baseline);
342             hgainmod->Fill(modId,gain);
343             hgain->Fill(gain);
344             if(!isgoodan){
345               hbadchan->SetBinContent(sideId+1,iddl+1, 1+hbadchan->GetBinContent(sideId+1,iddl+1));
346             }else{
347               nGoodAnodes++;
348             }
349           }
350           fclose(basFil);
351         }
352       }
353     }
354   }
355   hrawnoisemod->SetStats(0);
356   hcornoisemod->SetStats(0);
357   hbasemod->SetStats(0);
358   hgainmod->SetStats(0);
359   hbadchan->SetStats(0);
360   gStyle->SetPalette(1);
361   TString txtCountGood;
362   Float_t fracGood=100.*(Float_t)nGoodAnodes/(520.*260.);
363   txtCountGood.Form("Number of GoodAnodes = %d (%5.1f\%)\n",nGoodAnodes,fracGood);
364
365   TCanvas* c0=new TCanvas("c0","Bad Channels",800,900);
366   c0->SetRightMargin(0.14);
367   c0->SetBottomMargin(0.2);
368   hbadchan->Draw("colz"); 
369   hbadchan->GetXaxis()->SetTitle("Channel");
370   hbadchan->GetYaxis()->SetTitle("DDL");
371   hbadchan->GetXaxis()->SetTickLength(0);
372   hbadchan->GetYaxis()->SetTickLength(0);
373   c0->cd();
374   TLine** linv=new TLine*[12];
375   for(Int_t i=0;i<12;i++){
376     linv[i]=new TLine(i+0.75,-0.5,i+0.75,23.5);
377     linv[i]->SetLineColor(kGray+1);
378     linv[i]->Draw();
379   }
380   TLine** linh=new TLine*[24];
381   for(Int_t i=0;i<24;i++){
382     linh[i]=new TLine(-0.25,i+0.5,11.75,i+0.5);
383     linh[i]->SetLineColor(kGray+1);
384     linh[i]->Draw();
385   }
386   TLatex* tg=new TLatex(0.1,0.05,txtCountGood.Data());
387   tg->SetNDC();
388   tg->SetTextColor(4);
389   tg->SetTextSize(0.04);
390   tg->Draw();
391   c0->Update();
392  
393   TCanvas* c1=new TCanvas("c1","Baseline",1200,700);
394   c1->Divide(2,1);
395   c1->cd(1);
396   gPad->SetLeftMargin(0.14);
397   gPad->SetRightMargin(0.14);
398   hbase->Draw();
399   hbase->GetXaxis()->SetTitle("Baseline (ADC)");
400   c1->cd(2);
401   gPad->SetLeftMargin(0.14);
402   gPad->SetRightMargin(0.14);
403   hbasemod->Draw("colz");
404   hbasemod->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
405   hbasemod->GetYaxis()->SetTitle("Baseline (ADC)");
406   hbasemod->GetYaxis()->SetTitleOffset(1.35);
407
408   TCanvas* c2=new TCanvas("c2","Raw Noise",1200,700);
409   c2->Divide(2,1);
410   c2->cd(1);
411   gPad->SetLeftMargin(0.14);
412   gPad->SetRightMargin(0.14);
413   hrawnoise->Draw();
414   hrawnoise->GetXaxis()->SetTitle("Raw Noise (ADC)");
415   c2->cd(2);
416   gPad->SetLeftMargin(0.14);
417   gPad->SetRightMargin(0.14);
418   hrawnoisemod->Draw("colz");
419   hrawnoisemod->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
420   hrawnoisemod->GetYaxis()->SetTitle("Raw Noise (ADC)");
421   hrawnoisemod->GetYaxis()->SetTitleOffset(1.35);
422
423   TCanvas* c3=new TCanvas("c3","Corr Noise",1200,700);
424   c3->Divide(2,1);
425   c3->cd(1);
426   gPad->SetLeftMargin(0.14);
427   gPad->SetRightMargin(0.14);
428   hcornoise->Draw();
429   hcornoise->GetXaxis()->SetTitle("Raw Noise (ADC)");
430   c3->cd(2);
431   gPad->SetLeftMargin(0.14);
432   gPad->SetRightMargin(0.14);
433   hcornoisemod->Draw("colz");
434   hcornoisemod->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
435   hcornoisemod->GetYaxis()->SetTitle("Raw Noise (ADC)");
436   hcornoisemod->GetYaxis()->SetTitleOffset(1.35);
437
438   TCanvas* c4=new TCanvas("c4","Gain",1200,700);
439   c4->Divide(2,1);
440   c4->cd(1);
441   gPad->SetLeftMargin(0.14);
442   gPad->SetRightMargin(0.14);
443   hgain->Draw();
444   hgain->GetXaxis()->SetTitle("Gain (ADC/DAC)");
445   c4->cd(2);
446   gPad->SetLeftMargin(0.14);
447   gPad->SetRightMargin(0.14);
448   hgainmod->Draw("colz");
449   hgainmod->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
450   hgainmod->GetYaxis()->SetTitle("Gain (ADC/DAC)");
451   hgainmod->GetYaxis()->SetTitleOffset(1.35);
452
453   return;
454 }
455
456 void CheckCalibInterface::ShowInjector(){
457   ClearAll();
458
459   TH2F* hinjstatus=new TH2F("hinjstatus","Injector Status",24,-0.25,11.75,24,-0.5,23.5);
460   TGraph *vvsmod0=new TGraph(0);
461   TGraph *vvsmod1=new TGraph(0);
462   TGraph *poldegvsmod0=new TGraph(0); 
463   TGraph *poldegvsmod1=new TGraph(0); 
464   TGraph *anmaxvsmod0=new TGraph(0); 
465   TGraph *anmaxvsmod1=new TGraph(0); 
466   TGraph *dvcevsmod0=new TGraph(0);
467   TGraph *dvcevsmod1=new TGraph(0);
468   TGraph *dveevsmod0=new TGraph(0);
469   TGraph *dveevsmod1=new TGraph(0);
470   vvsmod0->SetTitle("Drift Speed vs. mod. number");
471   vvsmod1->SetTitle("Drift Speed vs. mod. number");
472   poldegvsmod0->SetTitle("Degree of poly fit vs. mod. number");
473   poldegvsmod1->SetTitle("Degree of poly fit vs. mod. number");
474   anmaxvsmod0->SetTitle("Anode with max. vdrift vs. mod. number");
475   anmaxvsmod1->SetTitle("Anode with max. vdrift vs. mod. number");
476   dvcevsmod0->SetTitle("Delta Vdrift 128-0 vs. mod. number");
477   dvcevsmod1->SetTitle("Delta Vdrift 128-0 vs. mod. number");
478   dveevsmod0->SetTitle("Delta Vdrift 256-0 vs. mod. number");
479   dveevsmod1->SetTitle("Delta Vdrift 256-0 vs. mod. number");
480
481   TF1* fPoly=new TF1("fPoly","pol3",0.,256.);
482
483   Int_t evNumb,polDeg; 
484   UInt_t timeStamp,statusInj;
485   Int_t retfscf;
486   TString inpFileName;
487   Float_t auxP;
488   Int_t iGoodInj=0;
489   for(Int_t iddl=0;iddl<24;iddl++){
490     for(Int_t imod=0;imod<12;imod++){
491       for(Int_t isid=0;isid<2;isid++){
492         inpFileName.Form("./calibFiles/SDDinj_ddl%02dc%02d_sid%d.data",iddl,imod,isid);
493         FILE* injFil = fopen(inpFileName.Data(),"read");
494         Int_t sideId=imod*2+isid;
495         Int_t modId=iddl*12+imod;
496         if (injFil == 0){ 
497           hinjstatus->SetBinContent(sideId+1,iddl+1,0);
498         }else{
499           Bool_t firstEvent=kTRUE;
500           retfscf=fscanf(injFil,"%d",&polDeg);
501           while (!feof(injFil)){
502             retfscf=fscanf(injFil,"%d %u ",&evNumb,&timeStamp);
503             if(evNumb==-99){
504               statusInj=timeStamp;
505               Int_t n7=(statusInj&(0x1F<<25))>>25;
506               Int_t n6=(statusInj&(0x1F<<20))>>20;
507               Int_t n5=(statusInj&(0x1F<<15))>>15;
508               Int_t n4=(statusInj&(0x1F<<5))>>10;
509               Int_t n3=(statusInj&(0x1F<<5))>>5;
510               Int_t n2=statusInj&0x1F;
511               Float_t aveStatus=(7.*n7+6.*n6+5.*n5+4.*n4+3.*n3+2.*n2)/32;
512               hinjstatus->SetBinContent(sideId+1,iddl+1,aveStatus);
513             }else{
514               if(feof(injFil)) break;
515               for(Int_t ic=0;ic<4;ic++){ 
516                 retfscf=fscanf(injFil,"%f ",&auxP);
517                 fPoly->SetParameter(ic,auxP);           
518               }
519               if(firstEvent && polDeg>0){
520                 firstEvent=kFALSE;
521                 ++iGoodInj;
522                 if(isid==0){
523                   vvsmod0->SetPoint(vvsmod0->GetN(),(Float_t)modId,fPoly->Eval(128));
524                   poldegvsmod0->SetPoint(poldegvsmod0->GetN(),(Float_t)modId,polDeg);
525                   anmaxvsmod0->SetPoint(anmaxvsmod0->GetN(),(Float_t)modId,fPoly->GetMaximumX(0.,256.));
526                   dvcevsmod0->SetPoint(dvcevsmod0->GetN(),(Float_t)modId,fPoly->Eval(128)-fPoly->Eval(0));
527                   dveevsmod0->SetPoint(dveevsmod0->GetN(),(Float_t)modId,fPoly->Eval(256)-fPoly->Eval(0));
528                 }else{
529                   vvsmod1->SetPoint(vvsmod1->GetN(),(Float_t)modId,fPoly->Eval(128));
530                   poldegvsmod1->SetPoint(poldegvsmod1->GetN(),(Float_t)modId,polDeg);
531                   anmaxvsmod1->SetPoint(anmaxvsmod1->GetN(),(Float_t)modId,fPoly->GetMaximumX(0.,256.));
532                   dvcevsmod1->SetPoint(dvcevsmod1->GetN(),(Float_t)modId,fPoly->Eval(128)-fPoly->Eval(0));
533                   dveevsmod1->SetPoint(dveevsmod1->GetN(),(Float_t)modId,fPoly->Eval(256)-fPoly->Eval(0));
534                 }
535               }
536             }
537           }     
538           fclose(injFil);
539         }
540       }
541     }
542   }
543   delete fPoly;
544
545   TString countmods;
546   countmods.Form("Number of half-modules with drift speed from injectors = %d",iGoodInj);
547   gStyle->SetPalette(51);
548   hinjstatus->SetStats(0);
549   hinjstatus->SetMinimum(-0.01);
550   hinjstatus->SetMaximum(7.);
551   TCanvas* c0=new TCanvas("c0","Injector status",800,900);
552   c0->SetRightMargin(0.14);
553   c0->SetBottomMargin(0.2);
554   hinjstatus->Draw("colz"); 
555   hinjstatus->GetXaxis()->SetTitle("Channel");
556   hinjstatus->GetYaxis()->SetTitle("DDL");
557   hinjstatus->GetXaxis()->SetTickLength(0);
558   hinjstatus->GetYaxis()->SetTickLength(0);
559   c0->cd();
560   TLine** linv=new TLine*[12];
561   for(Int_t i=0;i<12;i++){
562     linv[i]=new TLine(i+0.75,-0.5,i+0.75,23.5);
563     linv[i]->SetLineColor(kGray+1);
564     linv[i]->Draw();
565   }
566   TLine** linh=new TLine*[24];
567   for(Int_t i=0;i<24;i++){
568     linh[i]=new TLine(-0.25,i+0.5,11.75,i+0.5);
569     linh[i]->SetLineColor(kGray+1);
570     linh[i]->Draw();
571   }
572   TLatex* t3=new TLatex(0.1,0.05,countmods.Data());
573   t3->SetNDC();
574   t3->SetTextColor(4);
575   t3->SetTextSize(0.03);
576   t3->Draw();
577   c0->Update();
578
579   TLatex* tleft=new TLatex(0.2,0.82,"Side 0");
580   tleft->SetNDC();
581   tleft->SetTextColor(1);
582   TLatex* tright=new TLatex(0.2,0.75,"Side 1");
583   tright->SetNDC();
584   tright->SetTextColor(2);
585
586   TCanvas* c1;
587   c1=new TCanvas("c1","Vdrift vs. mod",1000,700);  
588   vvsmod0->SetMarkerStyle(20);
589   vvsmod0->Draw("AP");
590   vvsmod0->GetXaxis()->SetLimits(-1,290);
591   vvsmod0->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
592   vvsmod0->GetYaxis()->SetTitle("Vdrift (#mum/ns)");
593   vvsmod1->SetMarkerStyle(21);
594   vvsmod1->SetMarkerColor(2);
595   vvsmod1->Draw("SAMEP");
596   tleft->Draw();
597   tright->Draw();
598   TLine* ltop=new TLine(12*8-0.5,vvsmod0->GetYaxis()->GetXmin(),12*8-0.5,vvsmod0->GetYaxis()->GetXmax());
599   ltop->SetLineColor(4);
600   ltop->Draw();
601   TLatex* ttop=new TLatex(12*3.5,vvsmod0->GetYaxis()->GetXmin()+0.05,"TOP");
602   ttop->SetTextColor(4);
603   ttop->Draw();
604   TLine* lmed=new TLine(12*16-0.5,vvsmod0->GetYaxis()->GetXmin(),12*16-0.5,vvsmod0->GetYaxis()->GetXmax());
605   lmed->SetLineColor(4);
606   lmed->Draw();
607   TLatex* tmed=new TLatex(12*11.5,vvsmod0->GetYaxis()->GetXmin()+0.05,"MED");
608   tmed->SetTextColor(4);
609   tmed->Draw();
610   TLatex* tbot=new TLatex(12*19.5,vvsmod0->GetYaxis()->GetXmin()+0.05,"BOT");
611   tbot->SetTextColor(4);
612   tbot->Draw();
613
614
615   TCanvas* c2=new TCanvas("c2","Params vs. mod",900,900);
616   c2->Divide(2,2);
617   
618   c2->cd(1);
619   gPad->SetLeftMargin(0.14);
620   poldegvsmod0->SetMarkerStyle(20);
621   poldegvsmod0->Draw("AP");
622   poldegvsmod0->GetXaxis()->SetLimits(-1,290);
623   poldegvsmod0->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
624   poldegvsmod0->GetYaxis()->SetTitle("Degree of Polynomial fit");
625   poldegvsmod0->GetYaxis()->SetTitleOffset(1.4);
626   poldegvsmod1->SetMarkerStyle(21);
627   poldegvsmod1->SetMarkerColor(2);
628   poldegvsmod1->Draw("SAMEP");
629   tleft->Draw();
630   tright->Draw();
631   c2->cd(2);
632   gPad->SetLeftMargin(0.14);
633   anmaxvsmod0->SetMarkerStyle(20);
634   anmaxvsmod0->Draw("AP");
635   anmaxvsmod0->GetXaxis()->SetLimits(-1,290);
636   anmaxvsmod0->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
637   anmaxvsmod0->GetYaxis()->SetTitle("Anode with max. drift speed");
638   anmaxvsmod0->GetYaxis()->SetTitleOffset(1.4);
639   anmaxvsmod1->SetMarkerStyle(21);
640   anmaxvsmod1->SetMarkerColor(2);
641   anmaxvsmod1->Draw("SAMEP");
642   tleft->Draw();
643   tright->Draw();
644   c2->cd(3);
645   gPad->SetLeftMargin(0.14);
646   dvcevsmod0->SetMarkerStyle(20);
647   dvcevsmod0->Draw("AP");
648   dvcevsmod0->GetXaxis()->SetLimits(-1,290);
649   dvcevsmod0->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
650   dvcevsmod0->GetYaxis()->SetTitle("vdrift(anode128)-vdrift(anode0)");
651   dvcevsmod0->GetYaxis()->SetTitleOffset(1.4);
652   dvcevsmod1->SetMarkerStyle(21);
653   dvcevsmod1->SetMarkerColor(2);
654   dvcevsmod1->Draw("SAMEP");
655   tleft->Draw();
656   tright->Draw();
657   c2->cd(4);
658   gPad->SetLeftMargin(0.14);
659   dveevsmod0->SetMarkerStyle(20);
660   dveevsmod0->Draw("AP");
661   dveevsmod0->GetXaxis()->SetLimits(-1,290);
662   dveevsmod0->GetYaxis()->SetTitleOffset(1.4);
663   dveevsmod0->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
664   dveevsmod0->GetYaxis()->SetTitle("vdrift(anode256)-vdrift(anode0)");
665   dveevsmod1->SetMarkerStyle(21);
666   dveevsmod1->SetMarkerColor(2);
667   dveevsmod1->Draw("SAMEP");
668   tleft->Draw();
669   tright->Draw();
670
671   return;
672 }
673
674 void CheckCalibInterface::ShowDDL(Int_t iddl){
675   //
676   ClearAll();
677   TF1* fPoly=new TF1("fPoly","pol3",0.,256.);
678   TString inpFileName;
679   TGraph* gdrsp[24];
680   Int_t retfscf;
681   Int_t evNumb,polDeg; 
682   UInt_t timeStamp,statusInj;
683   Float_t auxP;
684   for(Int_t imod=0;imod<12;imod++){
685     for(Int_t isid=0;isid<2;isid++){
686       inpFileName.Form("./calibFiles/SDDinj_ddl%02dc%02d_sid%d.data",iddl,imod,isid);
687       FILE* injFil = fopen(inpFileName.Data(),"read");
688       Int_t sideId=imod*2+isid;
689       gdrsp[sideId]=new TGraph(0);
690       gdrsp[sideId]->SetTitle(Form("DDL %d  Mod %d  Sid %d\n",iddl,imod,isid));
691       Bool_t firstEvent=kTRUE;
692       if (injFil != 0){
693         retfscf=fscanf(injFil,"%d",&polDeg);
694         while (!feof(injFil)){
695           retfscf=fscanf(injFil,"%d %u ",&evNumb,&timeStamp);
696           if(evNumb==-99){
697             statusInj=timeStamp;
698           }else{
699             if(feof(injFil)) break;
700             for(Int_t ic=0;ic<4;ic++){ 
701               retfscf=fscanf(injFil,"%f ",&auxP);
702               fPoly->SetParameter(ic,auxP);
703             }     
704           }
705           if(firstEvent==kTRUE && polDeg>0){
706             firstEvent=kFALSE;
707             for(Int_t ian=0; ian<256; ian+=8){
708               gdrsp[sideId]->SetPoint(gdrsp[sideId]->GetN(),(Float_t)ian,fPoly->Eval(ian));
709             }
710           }
711         }
712         fclose(injFil); 
713       }
714       if(gdrsp[sideId]->GetN()==0) gdrsp[sideId]->SetPoint(0,128.,0.);
715     }
716   }
717   delete fPoly;
718
719   TCanvas* c0=new TCanvas("c0","Drift Speed Vs. Anode",1300,900);
720   c0->Divide(6,4);
721   for(Int_t i=0; i<24;i++){
722     c0->cd(i+1);
723     gdrsp[i]->Draw("APL");
724     gdrsp[i]->GetXaxis()->SetTitle("Anode");
725     gdrsp[i]->GetYaxis()->SetTitle("Drift Speed (#mum/ns)");
726     gdrsp[i]->GetYaxis()->SetTitleOffset(1.35);
727   }
728
729 }
730
731 void CheckCalibInterface::ShowSingleModule(Int_t iddl, Int_t ichan){
732   //
733   ClearAll();
734   TString inpFileName1,inpFileName2,inpFileName3;
735   Int_t retfscf;
736   Float_t baseline,rawnoise,cmn,corn,gain;
737   Int_t isgoodan,i,basmin,basoff;
738   Int_t th,tl;
739   Int_t iii,jjj,kkk;
740   Int_t evNumb,polDeg; 
741   UInt_t timeStamp,statusInj;
742   Float_t auxP;
743   TF1* fPoly=new TF1("fPoly","pol3",0.,256.);
744
745   TGraph* gbasel=new TGraph(0);
746   TGraph* gbaser=new TGraph(0);
747   TGraph* grawnl=new TGraph(0);
748   TGraph* grawnr=new TGraph(0);
749   TGraph* gcornl=new TGraph(0);
750   TGraph* gcornr=new TGraph(0);
751   TGraph* ggainl=new TGraph(0);
752   TGraph* ggainr=new TGraph(0);
753   TGraph* gstpdl=new TGraph(0);
754   TGraph* gstpdr=new TGraph(0);
755   TGraph* gstpul=new TGraph(0);
756   TGraph* gstpur=new TGraph(0);
757   TGraph* gdrspl=new TGraph(0);
758   TGraph* gdrspr=new TGraph(0);
759   gbasel->SetTitle("Baseline Left");
760   gbaser->SetTitle("Baseline Right");
761   grawnl->SetTitle("Noise Left");
762   grawnr->SetTitle("Noise Right");
763   gcornl->SetTitle("Noise Left");
764   gcornr->SetTitle("Noise Right");
765   ggainl->SetTitle("Gain Left");
766   ggainr->SetTitle("Gain Right");
767   gstpdl->SetTitle("Status Left");
768   gstpdr->SetTitle("Status Right");
769   gstpul->SetTitle("Status Left");
770   gstpur->SetTitle("Status Right");
771   gdrspl->SetTitle("Drift Speed Left");
772   gdrspr->SetTitle("Drift Speed Right");
773
774   for(Int_t isid=0;isid<2;isid++){
775     inpFileName1.Form("./calibFiles/SDDbase_step2_ddl%02dc%02d_sid%d.data",iddl,ichan,isid);
776     inpFileName2.Form("./calibFiles/SDDbase_ddl%02dc%02d_sid%d.data",iddl,ichan,isid);
777     inpFileName3.Form("./calibFiles/SDDinj_ddl%02dc%02d_sid%d.data",iddl,ichan,isid);
778   
779   
780     FILE* basFil = fopen(inpFileName1.Data(),"read");
781     if (basFil != 0){
782       retfscf=fscanf(basFil,"%d\n",&th);
783       retfscf=fscanf(basFil,"%d\n",&tl);
784       if(th==255 && tl==20) continue;
785       for(Int_t ian=0;ian<256;ian++){
786         retfscf=fscanf(basFil,"%d %d %f %d %d %f %f %f\n",&i,&isgoodan,&baseline,&basmin,&basoff,&rawnoise,&cmn,&corn);
787         if(isid==0){
788           gbasel->SetPoint(gbasel->GetN(),(Float_t)i,baseline);
789           grawnl->SetPoint(grawnl->GetN(),(Float_t)i,rawnoise);
790           gcornl->SetPoint(gcornl->GetN(),(Float_t)i,corn);
791           gstpdl->SetPoint(gstpdl->GetN(),(Float_t)i,(Float_t)isgoodan);
792         }else{
793           gbaser->SetPoint(gbaser->GetN(),(Float_t)i,baseline);
794           grawnr->SetPoint(grawnr->GetN(),(Float_t)i,rawnoise);
795           gcornr->SetPoint(gcornr->GetN(),(Float_t)i,corn);
796           gstpdr->SetPoint(gstpdr->GetN(),(Float_t)i,(Float_t)isgoodan);
797         }
798       }
799       fclose(basFil);
800     }
801         
802     FILE* pulFil = fopen(inpFileName2.Data(),"read");
803     if (pulFil != 0){
804       retfscf=fscanf(pulFil,"%d %d %d\n",&iii,&jjj,&kkk);
805       retfscf=fscanf(pulFil,"%d\n",&th);
806       retfscf=fscanf(pulFil,"%d\n",&tl);
807       if(th==255 && tl==20) continue;
808       for(Int_t ian=0;ian<256;ian++){
809         retfscf=fscanf(pulFil,"%d %d %f %d %d %f %f %f %f\n",&i,&isgoodan,&baseline,&basmin,&basoff,&rawnoise,&cmn,&corn,&gain);
810         if(isid==0){
811           ggainl->SetPoint(ggainl->GetN(),(Float_t)i,gain);
812           gstpul->SetPoint(gstpul->GetN(),(Float_t)i,(Float_t)isgoodan);
813         }else{
814           ggainr->SetPoint(ggainr->GetN(),(Float_t)i,gain);
815           gstpur->SetPoint(gstpur->GetN(),(Float_t)i,(Float_t)isgoodan);
816         }
817       }
818       fclose(pulFil);
819     }
820     
821     FILE* injFil = fopen(inpFileName3.Data(),"read");
822     Bool_t firstEvent=kTRUE;
823     if (injFil != 0){
824       retfscf=fscanf(injFil,"%d",&polDeg);
825       while (!feof(injFil)){
826         retfscf=fscanf(injFil,"%d %u ",&evNumb,&timeStamp);
827         if(evNumb==-99){
828           statusInj=timeStamp;
829         }else{
830           if(feof(injFil)) break;
831           for(Int_t ic=0;ic<4;ic++){ 
832             retfscf=fscanf(injFil,"%f ",&auxP);
833             fPoly->SetParameter(ic,auxP);
834           }       
835         }
836         if(firstEvent==kTRUE && polDeg>0){
837           firstEvent=kFALSE;
838           for(Int_t ian=0; ian<256; ian+=8){
839             if(isid==0) gdrspl->SetPoint(gdrspl->GetN(),(Float_t)ian,fPoly->Eval(ian));
840             else gdrspr->SetPoint(gdrspr->GetN(),(Float_t)ian,fPoly->Eval(ian));
841             
842           }
843         }
844       }
845     }
846   }
847
848
849   TCanvas * c0=new TCanvas("c0","Baselines",900,600);
850   c0->Divide(2,1);
851   c0->cd(1);
852   gPad->SetLeftMargin(0.14);
853   gbasel->SetMarkerStyle(7);
854   gbasel->Draw("AP");
855   gbasel->GetXaxis()->SetTitle("Anode");
856   gbasel->GetYaxis()->SetTitle("Baseline (ADC)");
857   gbasel->GetYaxis()->SetTitleOffset(1.35);
858   c0->cd(2);
859   gPad->SetLeftMargin(0.14);
860   gbaser->SetMarkerStyle(7);
861   gbaser->Draw("AP");
862   gbaser->GetXaxis()->SetTitle("Anode");
863   gbaser->GetYaxis()->SetTitle("Baseline (ADC)");
864   gbaser->GetYaxis()->SetTitleOffset(1.35);
865
866   TCanvas * c1=new TCanvas("c1","Noise",900,600);
867   c1->Divide(2,1);
868   c1->cd(1);
869   gPad->SetLeftMargin(0.14);
870   TLatex* t1=new TLatex(0.6,0.8,"Raw");
871   t1->SetNDC();
872   TLatex* t2=new TLatex(0.6,0.75,"Corrected");
873   t2->SetTextColor(2);
874   t2->SetNDC();
875   grawnl->SetMarkerStyle(7);
876   grawnl->Draw("AP");
877   grawnl->GetXaxis()->SetTitle("Anode");
878   grawnl->GetYaxis()->SetTitle("Noise (ADC)");
879   grawnl->GetYaxis()->SetTitleOffset(1.35);
880   gcornl->SetMarkerStyle(7);
881   gcornl->SetMarkerColor(2);
882   gcornl->Draw("SAMEP");
883   t1->Draw();
884   t2->Draw();
885   c1->cd(2);
886   gPad->SetLeftMargin(0.14);
887   grawnr->SetMarkerStyle(7);
888   grawnr->Draw("AP");
889   grawnr->GetXaxis()->SetTitle("Anode");
890   grawnr->GetYaxis()->SetTitle("Noise (ADC)");
891   grawnr->GetYaxis()->SetTitleOffset(1.35);
892   gcornr->SetMarkerStyle(7);
893   gcornr->SetMarkerColor(2);
894   gcornr->Draw("SAMEP");
895   t1->Draw();
896   t2->Draw();
897  
898   TCanvas * c2=new TCanvas("c2","Gain",900,600);
899   c2->Divide(2,1);
900   c2->cd(1);
901   gPad->SetLeftMargin(0.14);
902   ggainl->SetMarkerStyle(7);
903   ggainl->Draw("AP");
904   ggainl->GetXaxis()->SetTitle("Anode");
905   ggainl->GetYaxis()->SetTitle("Gain (ADC/DAC)");
906   ggainl->GetYaxis()->SetTitleOffset(1.35);
907   c2->cd(2);
908   gPad->SetLeftMargin(0.14);
909   ggainr->SetMarkerStyle(7);
910   ggainr->Draw("AP");
911   ggainr->GetXaxis()->SetTitle("Anode");
912   ggainr->GetYaxis()->SetTitle("Gain (ADC/DAC)");
913   ggainr->GetYaxis()->SetTitleOffset(1.35);
914
915   TCanvas * c3=new TCanvas("c3","Anode Status",900,600);
916   c3->Divide(2,1);
917   c3->cd(1);
918   gPad->SetLeftMargin(0.14);
919   TLatex* t3=new TLatex(0.6,0.85,"Pedestal");
920   t3->SetNDC();
921   TLatex* t4=new TLatex(0.6,0.8,"Pulser");
922   t4->SetTextColor(4);
923   t4->SetNDC();
924   gstpdl->SetMarkerStyle(7);
925   gstpdl->Draw("AP");
926   gstpdl->SetMinimum(-0.001);
927   gstpdl->SetMaximum(1.2);
928   gstpdl->GetXaxis()->SetTitle("Anode");
929   gstpdl->GetYaxis()->SetTitle("Status (1=OK, 0=BAD)");
930   gstpdl->GetYaxis()->SetTitleOffset(1.35);
931   gstpul->SetMarkerStyle(7);
932   gstpul->SetMarkerColor(4);
933   gstpul->Draw("SAMEP");
934   t3->Draw();
935   t4->Draw();
936   c3->cd(2);
937   gPad->SetLeftMargin(0.14);
938   gstpdr->SetMarkerStyle(7);
939   gstpdr->Draw("AP");
940   gstpdr->SetMinimum(-0.001);
941   gstpdr->SetMaximum(1.2);
942   gstpdr->GetXaxis()->SetTitle("Anode");
943   gstpdr->GetYaxis()->SetTitle("Status (1=OK, 0=BAD)");
944   gstpdr->GetYaxis()->SetTitleOffset(1.35);
945   gstpur->SetMarkerStyle(7);
946   gstpur->SetMarkerColor(4);
947   gstpur->Draw("SAMEP");
948   t3->Draw();
949   t4->Draw();
950
951   TCanvas * c4=new TCanvas("c4","Drift Speed",900,600);
952   c4->Divide(2,1);
953   c4->cd(1);
954   gPad->SetLeftMargin(0.14);
955   gdrspl->SetMarkerStyle(7);
956   gdrspl->Draw("APL");
957   gdrspl->GetXaxis()->SetTitle("Anode");
958   gdrspl->GetYaxis()->SetTitle("Drift Speed (#mum/ns)");
959   gdrspl->GetYaxis()->SetTitleOffset(1.35);
960   c4->cd(2);
961   gPad->SetLeftMargin(0.14);
962   gdrspr->SetMarkerStyle(7);
963   gdrspr->Draw("APL");
964   gdrspr->GetXaxis()->SetTitle("Anode");
965   gdrspr->GetYaxis()->SetTitle("Drift Speed (#mum/ns)");
966   gdrspr->GetYaxis()->SetTitleOffset(1.35);
967
968   return;
969 }
970
971 void CheckCalibInterface::DoSetlabel()
972 {
973    // Slot method connected to the ValueSet(Long_t) signal.
974    // It displays the value set in TGNumberEntry widget.
975   fNumDDL=fDDL->GetNumberEntry()->GetIntNumber();
976   fNumChannel=fChannel->GetNumberEntry()->GetIntNumber();
977   fShowOneMod->SetCommand(Form("CheckCalibInterface::ShowSingleModule(%d,%d)",
978                                fNumDDL,fNumChannel));
979   fShowDDL->SetCommand(Form("CheckCalibInterface::ShowDDL(%d)",
980                                fNumDDL));
981 }
982
983 void CheckCalibInterface::ClearAll(){
984
985   gROOT->Clear();
986   if(gROOT->FindObject("c0")) delete gROOT->FindObject("c0");
987   if(gROOT->FindObject("c1")) delete gROOT->FindObject("c1");
988   if(gROOT->FindObject("c2")) delete gROOT->FindObject("c2");
989   if(gROOT->FindObject("c3")) delete gROOT->FindObject("c3");
990   if(gROOT->FindObject("c4")) delete gROOT->FindObject("c4");
991   if(gROOT->FindObject("hbadchan")) delete gROOT->FindObject("hbadchan");
992   if(gROOT->FindObject("hrawnoisemod")) delete gROOT->FindObject("hrawnoisemod");
993   if(gROOT->FindObject("hrawnoise")) delete gROOT->FindObject("hrawnoise");
994   if(gROOT->FindObject("hcornoisemod")) delete gROOT->FindObject("hcornoisemod");
995   if(gROOT->FindObject("hcornoise")) delete gROOT->FindObject("hcornoise");
996   if(gROOT->FindObject("hbasemod")) delete gROOT->FindObject("hbasemod");
997   if(gROOT->FindObject("hbase")) delete gROOT->FindObject("hbase");
998   if(gROOT->FindObject("hgainmod")) delete gROOT->FindObject("hgainmod");
999   if(gROOT->FindObject("hgain")) delete gROOT->FindObject("hgainnoise");
1000
1001 }
1002
1003 void CheckCalibs()
1004 {
1005    new CheckCalibInterface(gClient->GetRoot(), 600, 400); 
1006 }