1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TApplication.h>
9 #include <TGNumberEntry.h>
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)
29 class CheckCalibInterface : public TGMainFrame {
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;
42 TGGroupFrame *fGframeALL;
43 TGGroupFrame *fGframeSING;
44 TGGroupFrame *fGframe1;
45 TGGroupFrame *fGframe2;
47 TGNumberEntry *fChannel;
52 CheckCalibInterface(const TGWindow *p, UInt_t w, UInt_t h);
53 virtual ~CheckCalibInterface();
56 static void CopyFiles();
57 static void ShowPedestal();
58 static void ShowPulser();
59 static void ShowInjector();
60 static void ShowSingleModule(Int_t iddl, Int_t ichan);
61 static void ClearAll();
63 ClassDef(CheckCalibInterface, 0)
66 CheckCalibInterface::CheckCalibInterface(const TGWindow *p, UInt_t w, UInt_t h)
67 : TGMainFrame(p, w, h)
70 fHor1 = new TGHorizontalFrame(this);
71 fHor2 = new TGHorizontalFrame(this);
72 fHor3 = new TGHorizontalFrame(this);
73 fHor4 = new TGHorizontalFrame(this);
75 TGLayoutHints *lh=new TGLayoutHints(kLHintsTop | kLHintsLeft, 5, 5, 5, 5);
76 TGLayoutHints *lhc=new TGLayoutHints(kLHintsCenterY | kLHintsExpandY, 5, 5, 5, 5);
77 TGLayoutHints *lh1=new TGLayoutHints(kLHintsTop | kLHintsExpandX, 5, 5, 5, 5);
79 fCopyFiles = new TGTextButton(fHor1, "&Copy Files", "CheckCalibInterface::CopyFiles()");
80 fHor1->AddFrame(fCopyFiles, lh1);
82 fGframeALL = new TGGroupFrame(fHor2,"All Modules",kHorizontalFrame);
83 fShowPedest = new TGTextButton(fGframeALL, "&Show Pedestal","CheckCalibInterface::ShowPedestal()");
84 fShowPulser = new TGTextButton(fGframeALL, "&Show Pulser","CheckCalibInterface::ShowPulser()");
85 fShowInject = new TGTextButton(fGframeALL, "&Show Injector","CheckCalibInterface::ShowInjector()");
86 fGframeALL->AddFrame(fShowPedest, lh1);
87 fGframeALL->AddFrame(fShowPulser, lh1);
88 fGframeALL->AddFrame(fShowInject, lh1);
89 fHor2->AddFrame(fGframeALL, lh1);
91 fGframeSING = new TGGroupFrame(fHor3,"Single Module",kHorizontalFrame);
92 fGframe1 = new TGGroupFrame(fGframeSING, "DDL");
93 fDDL = new TGNumberEntry(fGframe1, 0, 2,-1, TGNumberFormat::kNESInteger,
94 TGNumberFormat::kNEANonNegative,
95 TGNumberFormat::kNELLimitMinMax,
97 fGframe2 = new TGGroupFrame(fGframeSING, "Channel");
98 fChannel = new TGNumberEntry(fGframe2, 0, 2,-1, TGNumberFormat::kNESInteger,
99 TGNumberFormat::kNEANonNegative,
100 TGNumberFormat::kNELLimitMinMax,
102 fGframe1->AddFrame(fDDL,lh);
103 fGframe2->AddFrame(fChannel,lh);
104 fGframeSING->AddFrame(fGframe1,lh);
105 fGframeSING->AddFrame(fGframe2,lh);
107 fDDL->Connect("ValueSet(Long_t)", "CheckCalibInterface", this, "DoSetlabel()");
108 (fDDL->GetNumberEntry())->Connect("ReturnPressed()", "CheckCalibInterface", this, "DoSetlabel()");
110 fChannel->Connect("ValueSet(Long_t)", "CheckCalibInterface", this, "DoSetlabel()");
111 (fChannel->GetNumberEntry())->Connect("ReturnPressed()", "CheckCalibInterface", this, "DoSetlabel()");
112 fHor3->AddFrame(fGframeSING, lh1);
114 fShowOneMod = new TGTextButton(fGframeSING, "&Show Selected Module","CheckCalibInterface::ShowSingleModule(0,0)");
115 fGframeSING->AddFrame(fShowOneMod, lhc);
117 fExit = new TGTextButton(fHor4, "&Exit", "gApplication->Terminate(0)");
118 fHor4->AddFrame(fExit, lh1);
125 SetCleanup(kDeepCleanup);
126 SetWindowName("Main Control");
128 Resize(GetDefaultSize());
132 CheckCalibInterface::~CheckCalibInterface()
140 void CheckCalibInterface::CopyFiles(){
143 TString ldcName[6]={"aldaqpc083","aldaqpc084","aldaqpc041",
144 "aldaqpc082","aldaqpc085","aldaqpc086"};
145 gSystem->Exec("rm -rf calibFiles");
146 gSystem->Exec("mkdir calibFiles");
147 for(Int_t iLDC=0; iLDC<6; iLDC++){
148 Int_t firstDDL=iLDC*4;
149 Int_t lastDDL=iLDC*4+3;
150 command.Form("scp %s:/dateSite/ldc-SDD-%02d-%02d-0/work/SDDbase_step2_LDC.tar calibFiles/.",ldcName[iLDC].Data(),firstDDL,lastDDL);
151 gSystem->Exec(command.Data());
152 command.Form("scp %s:/dateSite/ldc-SDD-%02d-%02d-0/work/SDDbase_LDC.tar calibFiles/.",ldcName[iLDC].Data(),firstDDL,lastDDL);
153 gSystem->Exec(command.Data());
154 command.Form("scp %s:/dateSite/ldc-SDD-%02d-%02d-0/work/SDDinj_LDC.tar calibFiles/.",ldcName[iLDC].Data(),firstDDL,lastDDL);
155 gSystem->Exec(command.Data());
156 gSystem->Exec("cd calibFiles; tar xvf SDDbase_step2_LDC.tar; tar xvf SDDbase_LDC.tar; tar xvf SDDinj_LDC.tar; cd ..");
158 printf("-------------- DONE ---------------\n");
163 void CheckCalibInterface::ShowPedestal(){
166 TH2F* hbadchan=new TH2F("hbadchan","Number of bad channels",24,-0.25,11.75,24,-0.5,23.5);
167 TH2F* hrawnoisemod=new TH2F("hrawnoisemod","",288,-0.5,287.5,50,0.,10.);
168 TH1F* hrawnoise=new TH1F("hrawnoise","",100,0.,10.);
169 TH2F* hcornoisemod=new TH2F("hcornoisemod","",288,-0.5,287.5,50,0.,10.);
170 TH1F* hcornoise=new TH1F("hcornoise","",100,0.,10.);
171 TH2F* hbasemod=new TH2F("hbasemod","",288,-0.5,287.5,50,0.,150.);
172 TH1F* hbase=new TH1F("hbase","",100,0.,150.);
175 Float_t baseline,rawnoise,cmn,corn;
176 Int_t isgoodan,i,basmin,basoff;
179 for(Int_t iddl=0;iddl<24;iddl++){
180 for(Int_t imod=0;imod<12;imod++){
181 for(Int_t isid=0;isid<2;isid++){
182 inpFileName.Form("./calibFiles/SDDbase_step2_ddl%02dc%02d_sid%d.data",iddl,imod,isid);
183 FILE* basFil = fopen(inpFileName.Data(),"read");
184 Int_t sideId=imod*2+isid;
185 Int_t modId=iddl*12+imod;
186 if (basFil == 0) hbadchan->SetBinContent(sideId+1,iddl+1,256);
188 retfscf=fscanf(basFil,"%d\n",&th);
189 retfscf=fscanf(basFil,"%d\n",&tl);
190 if(th==255 && tl==20) continue;
191 for(Int_t ian=0;ian<256;ian++){
192 retfscf=fscanf(basFil,"%d %d %f %d %d %f %f %f\n",&i,&isgoodan,&baseline,&basmin,&basoff,&rawnoise,&cmn,&corn);
193 hrawnoisemod->Fill(modId,rawnoise);
194 hrawnoise->Fill(rawnoise);
195 hcornoisemod->Fill(modId,corn);
196 hcornoise->Fill(corn);
197 hbasemod->Fill(modId,baseline);
198 hbase->Fill(baseline);
200 hbadchan->SetBinContent(sideId+1,iddl+1, 1+hbadchan->GetBinContent(sideId+1,iddl+1));
210 hrawnoisemod->SetStats(0);
211 hcornoisemod->SetStats(0);
212 hbasemod->SetStats(0);
213 hbadchan->SetStats(0);
214 gStyle->SetPalette(1);
215 TString txtCountGood;
216 Float_t fracGood=100.*(Float_t)nGoodAnodes/(520.*260.);
217 txtCountGood.Form("Number of GoodAnodes = %d (%5.1f\%)\n",nGoodAnodes,fracGood);
219 TCanvas* c0=new TCanvas("c0","Bad Channels",800,900);
220 c0->SetRightMargin(0.14);
221 c0->SetBottomMargin(0.2);
222 hbadchan->Draw("colz");
223 hbadchan->GetXaxis()->SetTitle("Channel");
224 hbadchan->GetYaxis()->SetTitle("DDL");
225 hbadchan->GetXaxis()->SetTickLength(0);
226 hbadchan->GetYaxis()->SetTickLength(0);
228 TLine** linv=new TLine*[12];
229 for(Int_t i=0;i<12;i++){
230 linv[i]=new TLine(i+0.75,-0.5,i+0.75,23.5);
231 linv[i]->SetLineColor(kGray+1);
234 TLine** linh=new TLine*[24];
235 for(Int_t i=0;i<24;i++){
236 linh[i]=new TLine(-0.25,i+0.5,11.75,i+0.5);
237 linh[i]->SetLineColor(kGray+1);
240 TLatex* tg=new TLatex(0.1,0.05,txtCountGood.Data());
243 tg->SetTextSize(0.04);
247 TCanvas* c1=new TCanvas("c1","Baseline",1200,700);
250 gPad->SetLeftMargin(0.14);
251 gPad->SetRightMargin(0.14);
253 hbase->GetXaxis()->SetTitle("Baseline (ADC)");
255 gPad->SetLeftMargin(0.14);
256 gPad->SetRightMargin(0.14);
257 hbasemod->Draw("colz");
258 hbasemod->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
259 hbasemod->GetYaxis()->SetTitle("Baseline (ADC)");
260 hbasemod->GetYaxis()->SetTitleOffset(1.35);
262 TCanvas* c2=new TCanvas("c2","Raw Noise",1200,700);
265 gPad->SetLeftMargin(0.14);
266 gPad->SetRightMargin(0.14);
268 hrawnoise->GetXaxis()->SetTitle("Raw Noise (ADC)");
270 gPad->SetLeftMargin(0.14);
271 gPad->SetRightMargin(0.14);
272 hrawnoisemod->Draw("colz");
273 hrawnoisemod->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
274 hrawnoisemod->GetYaxis()->SetTitle("Raw Noise (ADC)");
275 hrawnoisemod->GetYaxis()->SetTitleOffset(1.35);
277 TCanvas* c3=new TCanvas("c3","Corr Noise",1200,700);
280 gPad->SetLeftMargin(0.14);
281 gPad->SetRightMargin(0.14);
283 hcornoise->GetXaxis()->SetTitle("Raw Noise (ADC)");
285 gPad->SetLeftMargin(0.14);
286 gPad->SetRightMargin(0.14);
287 hcornoisemod->Draw("colz");
288 hcornoisemod->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
289 hcornoisemod->GetYaxis()->SetTitle("Raw Noise (ADC)");
290 hcornoisemod->GetYaxis()->SetTitleOffset(1.35);
293 void CheckCalibInterface::ShowPulser(){
296 TH2F* hbadchan=new TH2F("hbadchan","Number of bad channels",24,-0.25,11.75,24,-0.5,23.5);
297 TH2F* hrawnoisemod=new TH2F("hrawnoisemod","",288,-0.5,287.5,50,0.,10.);
298 TH1F* hrawnoise=new TH1F("hrawnoise","",100,0.,10.);
299 TH2F* hcornoisemod=new TH2F("hcornoisemod","",288,-0.5,287.5,50,0.,10.);
300 TH1F* hcornoise=new TH1F("hcornoise","",100,0.,10.);
301 TH2F* hbasemod=new TH2F("hbasemod","",288,-0.5,287.5,50,0.,150.);
302 TH1F* hbase=new TH1F("hbase","",100,0.,150.);
303 TH2F* hgainmod=new TH2F("hgainmod","",288,-0.5,287.5,50,0.,5.);
304 TH1F* hgain=new TH1F("hgain","",100,0.,5.);
308 Float_t baseline,rawnoise,cmn,corn,gain;
309 Int_t isgoodan,i,basmin,basoff;
314 for(Int_t iddl=0;iddl<24;iddl++){
315 for(Int_t imod=0;imod<12;imod++){
316 for(Int_t isid=0;isid<2;isid++){
317 inpFileName.Form("./calibFiles/SDDbase_ddl%02dc%02d_sid%d.data",iddl,imod,isid);
318 FILE* basFil = fopen(inpFileName.Data(),"read");
319 Int_t sideId=imod*2+isid;
320 Int_t modId=iddl*12+imod;
321 if (basFil == 0) hbadchan->SetBinContent(sideId+1,iddl+1,256);
323 retfscf=fscanf(basFil,"%d %d %d\n",&iii,&jjj,&kkk);
325 retfscf=fscanf(basFil,"%d\n",&th);
326 retfscf=fscanf(basFil,"%d\n",&tl);
327 if(th==255 && tl==20) continue;
328 for(Int_t ian=0;ian<256;ian++){
329 retfscf=fscanf(basFil,"%d %d %f %d %d %f %f %f %f\n",&i,&isgoodan,&baseline,&basmin,&basoff,&rawnoise,&cmn,&corn,&gain);
330 hrawnoisemod->Fill(modId,rawnoise);
331 hrawnoise->Fill(rawnoise);
332 hcornoisemod->Fill(modId,corn);
333 hcornoise->Fill(corn);
334 hbasemod->Fill(modId,baseline);
335 hbase->Fill(baseline);
336 hgainmod->Fill(modId,gain);
339 hbadchan->SetBinContent(sideId+1,iddl+1, 1+hbadchan->GetBinContent(sideId+1,iddl+1));
349 hrawnoisemod->SetStats(0);
350 hcornoisemod->SetStats(0);
351 hbasemod->SetStats(0);
352 hgainmod->SetStats(0);
353 hbadchan->SetStats(0);
354 gStyle->SetPalette(1);
355 TString txtCountGood;
356 Float_t fracGood=100.*(Float_t)nGoodAnodes/(520.*260.);
357 txtCountGood.Form("Number of GoodAnodes = %d (%5.1f\%)\n",nGoodAnodes,fracGood);
359 TCanvas* c0=new TCanvas("c0","Bad Channels",800,900);
360 c0->SetRightMargin(0.14);
361 c0->SetBottomMargin(0.2);
362 hbadchan->Draw("colz");
363 hbadchan->GetXaxis()->SetTitle("Channel");
364 hbadchan->GetYaxis()->SetTitle("DDL");
365 hbadchan->GetXaxis()->SetTickLength(0);
366 hbadchan->GetYaxis()->SetTickLength(0);
368 TLine** linv=new TLine*[12];
369 for(Int_t i=0;i<12;i++){
370 linv[i]=new TLine(i+0.75,-0.5,i+0.75,23.5);
371 linv[i]->SetLineColor(kGray+1);
374 TLine** linh=new TLine*[24];
375 for(Int_t i=0;i<24;i++){
376 linh[i]=new TLine(-0.25,i+0.5,11.75,i+0.5);
377 linh[i]->SetLineColor(kGray+1);
380 TLatex* tg=new TLatex(0.1,0.05,txtCountGood.Data());
383 tg->SetTextSize(0.04);
387 TCanvas* c1=new TCanvas("c1","Baseline",1200,700);
390 gPad->SetLeftMargin(0.14);
391 gPad->SetRightMargin(0.14);
393 hbase->GetXaxis()->SetTitle("Baseline (ADC)");
395 gPad->SetLeftMargin(0.14);
396 gPad->SetRightMargin(0.14);
397 hbasemod->Draw("colz");
398 hbasemod->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
399 hbasemod->GetYaxis()->SetTitle("Baseline (ADC)");
400 hbasemod->GetYaxis()->SetTitleOffset(1.35);
402 TCanvas* c2=new TCanvas("c2","Raw Noise",1200,700);
405 gPad->SetLeftMargin(0.14);
406 gPad->SetRightMargin(0.14);
408 hrawnoise->GetXaxis()->SetTitle("Raw Noise (ADC)");
410 gPad->SetLeftMargin(0.14);
411 gPad->SetRightMargin(0.14);
412 hrawnoisemod->Draw("colz");
413 hrawnoisemod->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
414 hrawnoisemod->GetYaxis()->SetTitle("Raw Noise (ADC)");
415 hrawnoisemod->GetYaxis()->SetTitleOffset(1.35);
417 TCanvas* c3=new TCanvas("c3","Corr Noise",1200,700);
420 gPad->SetLeftMargin(0.14);
421 gPad->SetRightMargin(0.14);
423 hcornoise->GetXaxis()->SetTitle("Raw Noise (ADC)");
425 gPad->SetLeftMargin(0.14);
426 gPad->SetRightMargin(0.14);
427 hcornoisemod->Draw("colz");
428 hcornoisemod->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
429 hcornoisemod->GetYaxis()->SetTitle("Raw Noise (ADC)");
430 hcornoisemod->GetYaxis()->SetTitleOffset(1.35);
432 TCanvas* c4=new TCanvas("c4","Gain",1200,700);
435 gPad->SetLeftMargin(0.14);
436 gPad->SetRightMargin(0.14);
438 hgain->GetXaxis()->SetTitle("Gain (ADC/DAC)");
440 gPad->SetLeftMargin(0.14);
441 gPad->SetRightMargin(0.14);
442 hgainmod->Draw("colz");
443 hgainmod->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
444 hgainmod->GetYaxis()->SetTitle("Gain (ADC/DAC)");
445 hgainmod->GetYaxis()->SetTitleOffset(1.35);
450 void CheckCalibInterface::ShowInjector(){
453 TH2F* hinjstatus=new TH2F("hinjstatus","Injector Status",24,-0.25,11.75,24,-0.5,23.5);
454 TGraph *vvsmod0=new TGraph(0);
455 TGraph *vvsmod1=new TGraph(0);
456 TGraph *poldegvsmod0=new TGraph(0);
457 TGraph *poldegvsmod1=new TGraph(0);
458 TGraph *anmaxvsmod0=new TGraph(0);
459 TGraph *anmaxvsmod1=new TGraph(0);
460 TGraph *dvcevsmod0=new TGraph(0);
461 TGraph *dvcevsmod1=new TGraph(0);
462 TGraph *dveevsmod0=new TGraph(0);
463 TGraph *dveevsmod1=new TGraph(0);
464 vvsmod0->SetTitle("Drift Speed vs. mod. number");
465 vvsmod1->SetTitle("Drift Speed vs. mod. number");
466 poldegvsmod0->SetTitle("Degree of poly fit vs. mod. number");
467 poldegvsmod1->SetTitle("Degree of poly fit vs. mod. number");
468 anmaxvsmod0->SetTitle("Anode with max. vdrift vs. mod. number");
469 anmaxvsmod1->SetTitle("Anode with max. vdrift vs. mod. number");
470 dvcevsmod0->SetTitle("Delta Vdrift 128-0 vs. mod. number");
471 dvcevsmod1->SetTitle("Delta Vdrift 128-0 vs. mod. number");
472 dveevsmod0->SetTitle("Delta Vdrift 256-0 vs. mod. number");
473 dveevsmod1->SetTitle("Delta Vdrift 256-0 vs. mod. number");
475 TF1* fPoly=new TF1("fPoly","pol3",0.,256.);
478 UInt_t timeStamp,statusInj;
483 for(Int_t iddl=0;iddl<24;iddl++){
484 for(Int_t imod=0;imod<12;imod++){
485 for(Int_t isid=0;isid<2;isid++){
486 inpFileName.Form("./calibFiles/SDDinj_ddl%02dc%02d_sid%d.data",iddl,imod,isid);
487 FILE* injFil = fopen(inpFileName.Data(),"read");
488 Int_t sideId=imod*2+isid;
489 Int_t modId=iddl*12+imod;
491 hinjstatus->SetBinContent(sideId+1,iddl+1,0);
493 Bool_t firstEvent=kTRUE;
494 retfscf=fscanf(injFil,"%d",&polDeg);
495 while (!feof(injFil)){
496 retfscf=fscanf(injFil,"%d %u ",&evNumb,&timeStamp);
499 Int_t n7=(statusInj&(0x1F<<25))>>25;
500 Int_t n6=(statusInj&(0x1F<<20))>>20;
501 Int_t n5=(statusInj&(0x1F<<15))>>15;
502 Int_t n4=(statusInj&(0x1F<<5))>>10;
503 Int_t n3=(statusInj&(0x1F<<5))>>5;
504 Int_t n2=statusInj&0x1F;
505 Float_t aveStatus=(7.*n7+6.*n6+5.*n5+4.*n4+3.*n3+2.*n2)/32;
506 hinjstatus->SetBinContent(sideId+1,iddl+1,aveStatus);
508 if(feof(injFil)) break;
509 for(Int_t ic=0;ic<4;ic++){
510 retfscf=fscanf(injFil,"%f ",&auxP);
511 fPoly->SetParameter(ic,auxP);
513 if(firstEvent && polDeg>0){
517 vvsmod0->SetPoint(vvsmod0->GetN(),(Float_t)modId,fPoly->Eval(128));
518 poldegvsmod0->SetPoint(poldegvsmod0->GetN(),(Float_t)modId,polDeg);
519 anmaxvsmod0->SetPoint(anmaxvsmod0->GetN(),(Float_t)modId,fPoly->GetMaximumX(0.,256.));
520 dvcevsmod0->SetPoint(dvcevsmod0->GetN(),(Float_t)modId,fPoly->Eval(128)-fPoly->Eval(0));
521 dveevsmod0->SetPoint(dveevsmod0->GetN(),(Float_t)modId,fPoly->Eval(256)-fPoly->Eval(0));
523 vvsmod1->SetPoint(vvsmod1->GetN(),(Float_t)modId,fPoly->Eval(128));
524 poldegvsmod1->SetPoint(poldegvsmod1->GetN(),(Float_t)modId,polDeg);
525 anmaxvsmod1->SetPoint(anmaxvsmod1->GetN(),(Float_t)modId,fPoly->GetMaximumX(0.,256.));
526 dvcevsmod1->SetPoint(dvcevsmod1->GetN(),(Float_t)modId,fPoly->Eval(128)-fPoly->Eval(0));
527 dveevsmod1->SetPoint(dveevsmod1->GetN(),(Float_t)modId,fPoly->Eval(256)-fPoly->Eval(0));
540 countmods.Form("Number of half-modules with drift speed from injectors = %d",iGoodInj);
541 gStyle->SetPalette(59);
542 hinjstatus->SetStats(0);
543 hinjstatus->SetMinimum(-0.01);
544 hinjstatus->SetMaximum(7.);
545 TCanvas* c0=new TCanvas("c0","Injector status",800,900);
546 c0->SetRightMargin(0.14);
547 c0->SetBottomMargin(0.2);
548 hinjstatus->Draw("colz");
549 hinjstatus->GetXaxis()->SetTitle("Channel");
550 hinjstatus->GetYaxis()->SetTitle("DDL");
551 hinjstatus->GetXaxis()->SetTickLength(0);
552 hinjstatus->GetYaxis()->SetTickLength(0);
554 TLine** linv=new TLine*[12];
555 for(Int_t i=0;i<12;i++){
556 linv[i]=new TLine(i+0.75,-0.5,i+0.75,23.5);
557 linv[i]->SetLineColor(kGray+1);
560 TLine** linh=new TLine*[24];
561 for(Int_t i=0;i<24;i++){
562 linh[i]=new TLine(-0.25,i+0.5,11.75,i+0.5);
563 linh[i]->SetLineColor(kGray+1);
566 TLatex* t3=new TLatex(0.1,0.05,countmods.Data());
569 t3->SetTextSize(0.03);
573 TLatex* tleft=new TLatex(0.2,0.82,"Side 0");
575 tleft->SetTextColor(1);
576 TLatex* tright=new TLatex(0.2,0.75,"Side 1");
578 tright->SetTextColor(2);
581 c1=new TCanvas("c1","Vdrift vs. mod",1000,700);
582 vvsmod0->SetMarkerStyle(20);
584 vvsmod0->GetXaxis()->SetLimits(-1,290);
585 vvsmod0->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
586 vvsmod0->GetYaxis()->SetTitle("Vdrift (#mum/ns)");
587 vvsmod1->SetMarkerStyle(21);
588 vvsmod1->SetMarkerColor(2);
589 vvsmod1->Draw("SAMEP");
592 TLine* ltop=new TLine(12*8-0.5,vvsmod0->GetYaxis()->GetXmin(),12*8-0.5,vvsmod0->GetYaxis()->GetXmax());
593 ltop->SetLineColor(4);
595 TLatex* ttop=new TLatex(12*3.5,vvsmod0->GetYaxis()->GetXmin()+0.05,"TOP");
596 ttop->SetTextColor(4);
598 TLine* lmed=new TLine(12*16-0.5,vvsmod0->GetYaxis()->GetXmin(),12*16-0.5,vvsmod0->GetYaxis()->GetXmax());
599 lmed->SetLineColor(4);
601 TLatex* tmed=new TLatex(12*11.5,vvsmod0->GetYaxis()->GetXmin()+0.05,"MED");
602 tmed->SetTextColor(4);
604 TLatex* tbot=new TLatex(12*19.5,vvsmod0->GetYaxis()->GetXmin()+0.05,"BOT");
605 tbot->SetTextColor(4);
609 TCanvas* c2=new TCanvas("c2","Params vs. mod",900,900);
613 gPad->SetLeftMargin(0.14);
614 poldegvsmod0->SetMarkerStyle(20);
615 poldegvsmod0->Draw("AP");
616 poldegvsmod0->GetXaxis()->SetLimits(-1,290);
617 poldegvsmod0->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
618 poldegvsmod0->GetYaxis()->SetTitle("Degree of Polynomial fit");
619 poldegvsmod0->GetYaxis()->SetTitleOffset(1.4);
620 poldegvsmod1->SetMarkerStyle(21);
621 poldegvsmod1->SetMarkerColor(2);
622 poldegvsmod1->Draw("SAMEP");
626 gPad->SetLeftMargin(0.14);
627 anmaxvsmod0->SetMarkerStyle(20);
628 anmaxvsmod0->Draw("AP");
629 anmaxvsmod0->GetXaxis()->SetLimits(-1,290);
630 anmaxvsmod0->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
631 anmaxvsmod0->GetYaxis()->SetTitle("Anode with max. drift speed");
632 anmaxvsmod0->GetYaxis()->SetTitleOffset(1.4);
633 anmaxvsmod1->SetMarkerStyle(21);
634 anmaxvsmod1->SetMarkerColor(2);
635 anmaxvsmod1->Draw("SAMEP");
639 gPad->SetLeftMargin(0.14);
640 dvcevsmod0->SetMarkerStyle(20);
641 dvcevsmod0->Draw("AP");
642 dvcevsmod0->GetXaxis()->SetLimits(-1,290);
643 dvcevsmod0->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
644 dvcevsmod0->GetYaxis()->SetTitle("vdrift(anode128)-vdrift(anode0)");
645 dvcevsmod0->GetYaxis()->SetTitleOffset(1.4);
646 dvcevsmod1->SetMarkerStyle(21);
647 dvcevsmod1->SetMarkerColor(2);
648 dvcevsmod1->Draw("SAMEP");
652 gPad->SetLeftMargin(0.14);
653 dveevsmod0->SetMarkerStyle(20);
654 dveevsmod0->Draw("AP");
655 dveevsmod0->GetXaxis()->SetLimits(-1,290);
656 dveevsmod0->GetYaxis()->SetTitleOffset(1.4);
657 dveevsmod0->GetXaxis()->SetTitle("Module index (=DDL*12+Channel)");
658 dveevsmod0->GetYaxis()->SetTitle("vdrift(anode256)-vdrift(anode0)");
659 dveevsmod1->SetMarkerStyle(21);
660 dveevsmod1->SetMarkerColor(2);
661 dveevsmod1->Draw("SAMEP");
668 void CheckCalibInterface::ShowSingleModule(Int_t iddl, Int_t ichan){
671 TString inpFileName1,inpFileName2,inpFileName3;
673 Float_t baseline,rawnoise,cmn,corn,gain;
674 Int_t isgoodan,i,basmin,basoff;
678 UInt_t timeStamp,statusInj;
680 TF1* fPoly=new TF1("fPoly","pol3",0.,256.);
682 TGraph* gbasel=new TGraph(0);
683 TGraph* gbaser=new TGraph(0);
684 TGraph* grawnl=new TGraph(0);
685 TGraph* grawnr=new TGraph(0);
686 TGraph* gcornl=new TGraph(0);
687 TGraph* gcornr=new TGraph(0);
688 TGraph* ggainl=new TGraph(0);
689 TGraph* ggainr=new TGraph(0);
690 TGraph* gstpdl=new TGraph(0);
691 TGraph* gstpdr=new TGraph(0);
692 TGraph* gstpul=new TGraph(0);
693 TGraph* gstpur=new TGraph(0);
694 TGraph* gdrspl=new TGraph(0);
695 TGraph* gdrspr=new TGraph(0);
696 gbasel->SetTitle("Baseline Left");
697 gbaser->SetTitle("Baseline Right");
698 grawnl->SetTitle("Noise Left");
699 grawnr->SetTitle("Noise Right");
700 gcornl->SetTitle("Noise Left");
701 gcornr->SetTitle("Noise Right");
702 ggainl->SetTitle("Gain Left");
703 ggainr->SetTitle("Gain Right");
704 gstpdl->SetTitle("Status Left");
705 gstpdr->SetTitle("Status Right");
706 gstpul->SetTitle("Status Left");
707 gstpur->SetTitle("Status Right");
708 gdrspl->SetTitle("Drift Speed Left");
709 gdrspr->SetTitle("Drift Speed Right");
711 for(Int_t isid=0;isid<2;isid++){
712 inpFileName1.Form("./calibFiles/SDDbase_step2_ddl%02dc%02d_sid%d.data",iddl,ichan,isid);
713 inpFileName2.Form("./calibFiles/SDDbase_ddl%02dc%02d_sid%d.data",iddl,ichan,isid);
714 inpFileName3.Form("./calibFiles/SDDinj_ddl%02dc%02d_sid%d.data",iddl,ichan,isid);
717 FILE* basFil = fopen(inpFileName1.Data(),"read");
719 retfscf=fscanf(basFil,"%d\n",&th);
720 retfscf=fscanf(basFil,"%d\n",&tl);
721 if(th==255 && tl==20) continue;
722 for(Int_t ian=0;ian<256;ian++){
723 retfscf=fscanf(basFil,"%d %d %f %d %d %f %f %f\n",&i,&isgoodan,&baseline,&basmin,&basoff,&rawnoise,&cmn,&corn);
725 gbasel->SetPoint(gbasel->GetN(),(Float_t)i,baseline);
726 grawnl->SetPoint(grawnl->GetN(),(Float_t)i,rawnoise);
727 gcornl->SetPoint(gcornl->GetN(),(Float_t)i,corn);
728 gstpdl->SetPoint(gstpdl->GetN(),(Float_t)i,(Float_t)isgoodan);
730 gbaser->SetPoint(gbaser->GetN(),(Float_t)i,baseline);
731 grawnr->SetPoint(grawnr->GetN(),(Float_t)i,rawnoise);
732 gcornr->SetPoint(gcornr->GetN(),(Float_t)i,corn);
733 gstpdr->SetPoint(gstpdr->GetN(),(Float_t)i,(Float_t)isgoodan);
739 FILE* pulFil = fopen(inpFileName2.Data(),"read");
741 retfscf=fscanf(pulFil,"%d %d %d\n",&iii,&jjj,&kkk);
742 retfscf=fscanf(pulFil,"%d\n",&th);
743 retfscf=fscanf(pulFil,"%d\n",&tl);
744 if(th==255 && tl==20) continue;
745 for(Int_t ian=0;ian<256;ian++){
746 retfscf=fscanf(pulFil,"%d %d %f %d %d %f %f %f %f\n",&i,&isgoodan,&baseline,&basmin,&basoff,&rawnoise,&cmn,&corn,&gain);
748 ggainl->SetPoint(ggainl->GetN(),(Float_t)i,gain);
749 gstpul->SetPoint(gstpul->GetN(),(Float_t)i,(Float_t)isgoodan);
751 ggainr->SetPoint(ggainr->GetN(),(Float_t)i,gain);
752 gstpur->SetPoint(gstpur->GetN(),(Float_t)i,(Float_t)isgoodan);
758 FILE* injFil = fopen(inpFileName3.Data(),"read");
759 Bool_t firstEvent=kTRUE;
761 retfscf=fscanf(injFil,"%d",&polDeg);
762 while (!feof(injFil)){
763 retfscf=fscanf(injFil,"%d %u ",&evNumb,&timeStamp);
767 if(feof(injFil)) break;
768 for(Int_t ic=0;ic<4;ic++){
769 retfscf=fscanf(injFil,"%f ",&auxP);
770 fPoly->SetParameter(ic,auxP);
773 if(firstEvent==kTRUE && polDeg>0){
775 for(Int_t ian=0; ian<256; ian+=8){
776 if(isid==0) gdrspl->SetPoint(gdrspl->GetN(),(Float_t)ian,fPoly->Eval(ian));
777 else gdrspr->SetPoint(gdrspr->GetN(),(Float_t)ian,fPoly->Eval(ian));
786 TCanvas * c0=new TCanvas("c0","Baselines",900,600);
789 gPad->SetLeftMargin(0.14);
790 gbasel->SetMarkerStyle(7);
792 gbasel->GetXaxis()->SetTitle("Anode");
793 gbasel->GetYaxis()->SetTitle("Baseline (ADC)");
794 gbasel->GetYaxis()->SetTitleOffset(1.35);
796 gPad->SetLeftMargin(0.14);
797 gbaser->SetMarkerStyle(7);
799 gbaser->GetXaxis()->SetTitle("Anode");
800 gbaser->GetYaxis()->SetTitle("Baseline (ADC)");
801 gbaser->GetYaxis()->SetTitleOffset(1.35);
803 TCanvas * c1=new TCanvas("c1","Noise",900,600);
806 gPad->SetLeftMargin(0.14);
807 TLatex* t1=new TLatex(0.6,0.8,"Raw");
809 TLatex* t2=new TLatex(0.6,0.75,"Corrected");
812 grawnl->SetMarkerStyle(7);
814 grawnl->GetXaxis()->SetTitle("Anode");
815 grawnl->GetYaxis()->SetTitle("Noise (ADC)");
816 grawnl->GetYaxis()->SetTitleOffset(1.35);
817 gcornl->SetMarkerStyle(7);
818 gcornl->SetMarkerColor(2);
819 gcornl->Draw("SAMEP");
823 gPad->SetLeftMargin(0.14);
824 grawnr->SetMarkerStyle(7);
826 grawnr->GetXaxis()->SetTitle("Anode");
827 grawnr->GetYaxis()->SetTitle("Noise (ADC)");
828 grawnr->GetYaxis()->SetTitleOffset(1.35);
829 gcornr->SetMarkerStyle(7);
830 gcornr->SetMarkerColor(2);
831 gcornr->Draw("SAMEP");
835 TCanvas * c2=new TCanvas("c2","Gain",900,600);
838 gPad->SetLeftMargin(0.14);
839 ggainl->SetMarkerStyle(7);
841 ggainl->GetXaxis()->SetTitle("Anode");
842 ggainl->GetYaxis()->SetTitle("Gain (ADC/DAC)");
843 ggainl->GetYaxis()->SetTitleOffset(1.35);
845 gPad->SetLeftMargin(0.14);
846 ggainr->SetMarkerStyle(7);
848 ggainr->GetXaxis()->SetTitle("Anode");
849 ggainr->GetYaxis()->SetTitle("Gain (ADC/DAC)");
850 ggainr->GetYaxis()->SetTitleOffset(1.35);
852 TCanvas * c3=new TCanvas("c3","Anode Status",900,600);
855 gPad->SetLeftMargin(0.14);
856 TLatex* t3=new TLatex(0.6,0.85,"Pedestal");
858 TLatex* t4=new TLatex(0.6,0.8,"Pulser");
861 gstpdl->SetMarkerStyle(7);
863 gstpdl->SetMinimum(-0.001);
864 gstpdl->SetMaximum(1.2);
865 gstpdl->GetXaxis()->SetTitle("Anode");
866 gstpdl->GetYaxis()->SetTitle("Status (1=OK, 0=BAD)");
867 gstpdl->GetYaxis()->SetTitleOffset(1.35);
868 gstpul->SetMarkerStyle(7);
869 gstpul->SetMarkerColor(4);
870 gstpul->Draw("SAMEP");
874 gPad->SetLeftMargin(0.14);
875 gstpdr->SetMarkerStyle(7);
877 gstpdr->SetMinimum(-0.001);
878 gstpdr->SetMaximum(1.2);
879 gstpdr->GetXaxis()->SetTitle("Anode");
880 gstpdr->GetYaxis()->SetTitle("Status (1=OK, 0=BAD)");
881 gstpdr->GetYaxis()->SetTitleOffset(1.35);
882 gstpur->SetMarkerStyle(7);
883 gstpur->SetMarkerColor(4);
884 gstpur->Draw("SAMEP");
888 TCanvas * c4=new TCanvas("c4","Drift Speed",900,600);
891 gPad->SetLeftMargin(0.14);
892 gdrspl->SetMarkerStyle(7);
894 gdrspl->GetXaxis()->SetTitle("Anode");
895 gdrspl->GetYaxis()->SetTitle("Drift Speed (#mum/ns)");
896 gdrspl->GetYaxis()->SetTitleOffset(1.35);
898 gPad->SetLeftMargin(0.14);
899 gdrspr->SetMarkerStyle(7);
901 gdrspr->GetXaxis()->SetTitle("Anode");
902 gdrspr->GetYaxis()->SetTitle("Drift Speed (#mum/ns)");
903 gdrspr->GetYaxis()->SetTitleOffset(1.35);
908 void CheckCalibInterface::DoSetlabel()
910 // Slot method connected to the ValueSet(Long_t) signal.
911 // It displays the value set in TGNumberEntry widget.
912 fNumDDL=fDDL->GetNumberEntry()->GetIntNumber();
913 fNumChannel=fChannel->GetNumberEntry()->GetIntNumber();
914 fShowOneMod->SetCommand(Form("CheckCalibInterface::ShowSingleModule(%d,%d)",
915 fNumDDL,fNumChannel));
918 void CheckCalibInterface::ClearAll(){
921 if(gROOT->FindObject("c0")) delete gROOT->FindObject("c0");
922 if(gROOT->FindObject("c1")) delete gROOT->FindObject("c1");
923 if(gROOT->FindObject("c2")) delete gROOT->FindObject("c2");
924 if(gROOT->FindObject("c3")) delete gROOT->FindObject("c3");
925 if(gROOT->FindObject("c4")) delete gROOT->FindObject("c4");
926 if(gROOT->FindObject("hbadchan")) delete gROOT->FindObject("hbadchan");
927 if(gROOT->FindObject("hrawnoisemod")) delete gROOT->FindObject("hrawnoisemod");
928 if(gROOT->FindObject("hrawnoise")) delete gROOT->FindObject("hrawnoise");
929 if(gROOT->FindObject("hcornoisemod")) delete gROOT->FindObject("hcornoisemod");
930 if(gROOT->FindObject("hcornoise")) delete gROOT->FindObject("hcornoise");
931 if(gROOT->FindObject("hbasemod")) delete gROOT->FindObject("hbasemod");
932 if(gROOT->FindObject("hbase")) delete gROOT->FindObject("hbase");
933 if(gROOT->FindObject("hgainmod")) delete gROOT->FindObject("hgainmod");
934 if(gROOT->FindObject("hgain")) delete gROOT->FindObject("hgainnoise");
940 new CheckCalibInterface(gClient->GetRoot(), 400, 400);