/[PAMELA software]/quicklook/tracker/flight/macros/FTrkCalibQLook_EXPERT.cxx
ViewVC logotype

Contents of /quicklook/tracker/flight/macros/FTrkCalibQLook_EXPERT.cxx

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.9 - (show annotations) (download)
Fri Aug 11 10:24:09 2006 UTC (18 years, 3 months ago) by pam-fi
Branch: MAIN
CVS Tags: v1r08
Changes since 1.8: +147 -31 lines
release v1r08

1 /**
2 * FTrkCalibQLookExpert.cxx
3 *
4 * autor: D.Fedele
5 * version v1r08
6 * Parameters:
7 * file - the data file to analyze
8 * step - select =1 in order to analyze one event at time
9 * fromevent - first event to analyze
10 * toevent - last event to analyze
11 * outdir - total path of output file
12 * outfile - extension of output file (pdf,ps,gif,jpg)
13 *
14 */
15 //
16 #include <iostream>
17 //
18 #include <TPaveText.h>
19 #include <TLatex.h>
20 #include <TCanvas.h>
21 #include <TFile.h>
22 #include <TTree.h>
23 #include <TGraph.h>
24 #include <TStyle.h>
25 #include <TString.h>
26 //
27 #include <PscuHeader.h>
28 #include <EventHeader.h>
29 #include <CalibTrk1Event.h>
30 #include <CalibTrk2Event.h>
31 //
32
33 typedef struct caltrk_def{
34 Int_t good0[2];
35 Int_t daqmode[12];
36 Int_t dspnum[12];
37 Int_t calibnum[12];
38 Int_t ncalev[12];
39 Int_t calfl[12];
40 Int_t ped1[12];
41 Int_t ped2[12];
42 Int_t ped3[12];
43 Int_t sig1[12];
44 Int_t sig2[12];
45 Int_t sig3[12];
46 Int_t nbad1[12];
47 Int_t nbad2[12];
48 Int_t nbad3[12];
49 Float_t dspped[12][3702];
50 Float_t dspsig[12][3702];
51 Float_t dspbad[12][3702];
52 Int_t crc_c[12][3];
53 Int_t crc_hc[12];
54 };
55
56
57 void FTrkCalibQLook_EXPERT(TString file,Int_t step,Int_t fromevent, Int_t toevent,TString outdir,TString outfile)
58 {
59 //
60 // obtain information about the data file and select the output dir
61 Int_t dwpos = file.Last('/');
62 Int_t dwpos1 = file.Last('.');
63 TString base,ffile ;
64 ffile=file(dwpos+1,dwpos1-(dwpos+1));
65 if(dwpos>0) base=file(0,dwpos);
66
67 TString out;
68 if(outdir.Length()==0){
69 out = base;
70 }else{
71 out = outdir;
72 }
73 if(out.Last('/')+1<out.Length()) out+="/";
74
75 //
76 // inizialise the variables and open the file
77 Int_t nevents=0;
78 Int_t minevent = 0;
79 Int_t maxevent = 0;
80 ULong64_t FOBT[2];
81
82 FOBT[0]=0;
83 FOBT[1]=0;
84
85 struct caltrk_def ctrk;
86 ctrk.good0[0]=0;
87 ctrk.good0[1]=0;
88 for(Int_t i=0;i<12;i++){
89 ctrk.daqmode[i]=0;
90 ctrk.dspnum[i]=0;
91 ctrk.calibnum[i]=0;
92 ctrk.ncalev[i]=0;
93 ctrk.calfl[i]=0;
94 ctrk.ped1[i]=0;
95 ctrk.ped2[i]=0;
96 ctrk.ped3[i]=0;
97 ctrk.sig1[i]=0;
98 ctrk.sig2[i]=0;
99 ctrk.sig3[i]=0;
100 ctrk.nbad1[i]=0;
101 ctrk.nbad2[i]=0;
102 ctrk.nbad3[i]=0;
103 ctrk.crc_c[i][0]=0;
104 ctrk.crc_c[i][1]=0;
105 ctrk.crc_c[i][2]=0;
106 ctrk.crc_hc[i]=0;
107 for(Int_t iii=0;iii<3072;iii++){
108 ctrk.dspped[i][iii]=0;
109 ctrk.dspsig[i][iii]=0;
110 ctrk.dspbad[i][iii]=0;
111 };
112
113 };
114
115 TFile *calibFile = new TFile(file);
116 if ( !calibFile ){
117 printf("No data file, exiting...\n");
118 return;
119 };
120 TTree *otr1,*otr2,*hotr,*totr;
121 pamela::EventHeader *eh1=0,*eh2=0,*eh4=0,*eh3=0;
122 pamela::PscuHeader *ph=0;
123 pamela::CalibTrk1Event *trk1 = 0;
124 pamela::CalibTrk2Event *trk2 = 0;
125 pamela::EventCounter *cod=0;
126
127 pamela::PacketType *pctp=0;
128
129 hotr = (TTree*)calibFile->Get("CalibHeader");
130 hotr->SetBranchAddress("Header", &eh4);
131
132 totr = (TTree*)calibFile->Get("CalibTrailer");
133 totr->SetBranchAddress("Header", &eh3);
134
135 otr1 = (TTree*)calibFile->Get("CalibTrk1");
136 otr1->SetBranchAddress("CalibTrk1", &trk1);
137 otr1->SetBranchAddress("Header",&eh1);
138 otr2 = (TTree*)calibFile->Get("CalibTrk2");
139 otr2->SetBranchAddress("CalibTrk2", &trk2);
140 otr2->SetBranchAddress("Header",&eh2);
141
142 if(otr1->GetEntries()==otr2->GetEntries())
143 nevents = otr1->GetEntries();
144 else{
145 printf("WARNING: CalibTrk1 entries is different from CalibTrk2 entries");
146 return;}
147 if (nevents<=0) {
148 calibFile->Close();
149 printf("No calibration packets found, exiting...\n");
150 return;
151 };
152
153 printf("Number of calibration packets: %i\n",nevents);
154
155 if ( fromevent > toevent && toevent ){
156 printf("It must be fromevent < toevent \n");
157 return;
158 };
159
160 if ( fromevent > nevents || fromevent < 0 ) {
161 printf("You can choose fromevent between 0 (all) and %i \n",nevents);
162 return;
163 };
164
165 if ( toevent > nevents || toevent < 0 ) {
166 printf("You can choose toevent between 0 (all) and %i \n",nevents);
167 return;
168 };
169 if ( fromevent == 0 ) {
170 minevent = 0;
171 maxevent = nevents;
172 } else {
173 minevent = fromevent - 1;
174 if ( toevent > 0 ){
175 maxevent = toevent;
176 } else if (toevent > nevents) {
177 maxevent = nevents;
178 } else {
179 maxevent = fromevent;
180 };
181 };
182
183 //**********************************************************************
184 //
185 // LOOP OVER EVENTS
186 //
187 //**********************************************************************
188
189 Int_t hcevent=hotr->GetEntries();
190 ULong64_t HOBT[hcevent], TOBT[hcevent];
191 for (Int_t i = 0; i < hcevent; i++){
192 totr->GetEntry(i);
193 hotr->GetEntry(i);
194 ph = eh4->GetPscuHeader();
195 HOBT[i]= ph->GetOrbitalTime();
196 ph = eh3->GetPscuHeader();
197 TOBT[i]= ph->GetOrbitalTime();
198 }
199
200 printf("\n Scan of calibration packets from %i to %i ... \n\n",minevent+1,maxevent);
201 for (Int_t i = minevent; i < maxevent; i++){
202
203 otr1->GetEntry(i);
204 otr2->GetEntry(i);
205
206 ctrk.good0[0]=trk1->good0;
207 ctrk.good0[1]=trk2->good0;
208 for (Int_t m = 0; m < 6; m++){
209 ph = eh1->GetPscuHeader();
210 cod = eh1->GetCounter();
211 FOBT[0]= ph->GetOrbitalTime();
212 ctrk.daqmode[trk1->DSPnumber[m]-1]=trk1->DAQmode[m];
213 ctrk.dspnum[trk1->DSPnumber[m]-1]=trk1->DSPnumber[m];
214 ctrk.calibnum[trk1->DSPnumber[m]-1]=trk1->calibnumber[m];
215 ctrk.ncalev[trk1->DSPnumber[m]-1]=trk1->ncalib_event[m];
216 ctrk.ped1[trk1->DSPnumber[m]-1]=trk1->ped_l1[m];
217 ctrk.ped2[trk1->DSPnumber[m]-1]=trk1->ped_l2[m];
218 ctrk.ped3[trk1->DSPnumber[m]-1]=trk1->ped_l3[m];
219 ctrk.sig1[trk1->DSPnumber[m]-1]=trk1->sig_l1[m];
220 ctrk.sig2[trk1->DSPnumber[m]-1]=trk1->sig_l2[m];
221 ctrk.sig3[trk1->DSPnumber[m]-1]=trk1->sig_l3[m];
222 ctrk.nbad1[trk1->DSPnumber[m]-1]=trk1->nbad_l1[m];
223 ctrk.nbad2[trk1->DSPnumber[m]-1]=trk1->nbad_l2[m];
224 ctrk.nbad3[trk1->DSPnumber[m]-1]=trk1->nbad_l3[m];
225 ctrk.calfl[trk1->DSPnumber[m]-1]=trk1->cal_flag[m];
226 ctrk.crc_c[trk1->DSPnumber[m]-1][0]=trk1->crc_cal[m][0];
227 ctrk.crc_c[trk1->DSPnumber[m]-1][1]=trk1->crc_cal[m][1];
228 ctrk.crc_c[trk1->DSPnumber[m]-1][2]=trk1->crc_cal[m][2];
229 ctrk.crc_hc[trk1->DSPnumber[m]-1]=trk1->crc_hcal[m];
230 for (Int_t j = 0; j < 3072; j++){
231 ctrk.dspped[trk1->DSPnumber[m]-1][j]=trk1->DSPped_par[m][j];
232 ctrk.dspsig[trk1->DSPnumber[m]-1][j]=trk1->DSPsig_par[m][j];
233 ctrk.dspbad[trk1->DSPnumber[m]-1][j]=trk1->DSPbad_par[m][j];
234 };
235 ph = eh2->GetPscuHeader();
236 FOBT[1]= ph->GetOrbitalTime();
237 ctrk.daqmode[trk2->DSPnumber[m]-1]=trk2->DAQmode[m];
238 ctrk.dspnum[trk2->DSPnumber[m]-1]=trk2->DSPnumber[m];
239 ctrk.calibnum[trk2->DSPnumber[m]-1]=trk2->calibnumber[m];
240 ctrk.ncalev[trk2->DSPnumber[m]-1]=trk2->ncalib_event[m];
241 ctrk.ped1[trk2->DSPnumber[m]-1]=trk2->ped_l1[m];
242 ctrk.ped2[trk2->DSPnumber[m]-1]=trk2->ped_l2[m];
243 ctrk.ped3[trk2->DSPnumber[m]-1]=trk2->ped_l3[m];
244 ctrk.sig1[trk2->DSPnumber[m]-1]=trk2->sig_l1[m];
245 ctrk.sig2[trk2->DSPnumber[m]-1]=trk2->sig_l2[m];
246 ctrk.sig3[trk2->DSPnumber[m]-1]=trk2->sig_l3[m];
247 ctrk.nbad1[trk2->DSPnumber[m]-1]=trk2->nbad_l1[m];
248 ctrk.nbad2[trk2->DSPnumber[m]-1]=trk2->nbad_l2[m];
249 ctrk.nbad3[trk2->DSPnumber[m]-1]=trk2->nbad_l3[m];
250 ctrk.calfl[trk2->DSPnumber[m]-1]=trk2->cal_flag[m];
251 ctrk.crc_c[trk1->DSPnumber[m]-1][0]=trk2->crc_cal[m][0];
252 ctrk.crc_c[trk1->DSPnumber[m]-1][1]=trk2->crc_cal[m][1];
253 ctrk.crc_c[trk1->DSPnumber[m]-1][2]=trk2->crc_cal[m][2];
254 ctrk.crc_hc[trk1->DSPnumber[m]-1]=trk2->crc_hcal[m];
255 for (Int_t j = 0; j < 3072; j++){
256 ctrk.dspped[trk2->DSPnumber[m]-1][j]=trk2->DSPped_par[m][j];
257 ctrk.dspsig[trk2->DSPnumber[m]-1][j]=trk2->DSPsig_par[m][j];
258 ctrk.dspbad[trk2->DSPnumber[m]-1][j]=trk2->DSPbad_par[m][j];
259 }
260 }
261
262
263 //
264 // other variables definitions
265 Int_t risposta=0;
266 stringstream fromfile;
267
268 fromfile<<"FTrkCalibQLook_EXPERT File: "<<ffile<<" -- CalibHeader OBT= "<<HOBT[(cod->Get(pctp->CalibHeader))-1]<<" -- Calib pkt OBT= "<<FOBT[0]<<" -- CalibTrailer OBT= "<<TOBT[(cod->Get(pctp->CalibHeader))-1]<<" --";
269
270 gStyle->SetLabelSize(0.07,"x");
271 gStyle->SetLabelSize(0.07,"y");
272 gStyle->SetTitleFillColor(10);
273 gStyle->SetTitleFontSize(0.08);
274 gStyle->SetTitleOffset(0.8,"y");
275 gStyle->SetTitleOffset(0.9,"x");
276 gStyle->SetTitleSize(0.06,"y");
277 gStyle->SetTitleSize(0.06,"x");
278 gStyle->SetOptStat(101110);
279 gStyle->SetStatX(0.9);
280 gStyle->SetStatW(0.4);
281 gStyle->SetStatColor(10);
282 gStyle->SetStatFontSize(0.1);
283
284 //
285 // draw display area
286
287 TLatex *tzz=new TLatex();
288 tzz->SetTextFont(32);
289 tzz->SetTextColor(1);
290 tzz->SetTextAlign(12);
291 tzz->SetTextSize(0.019);
292 Int_t canvasx=1200;
293 Int_t canvasy=900;
294 TCanvas *c1 = new TCanvas("c1","FTrkCalibQLook_EXPERT_ped",canvasx,canvasy);
295 c1->SetFillColor(10);
296 tzz->DrawLatex(.01,0.98,fromfile.str().c_str());
297 tzz->DrawLatex(.90,0.98,"PEDESTAL");
298 TCanvas *c2 = new TCanvas("c2","FTrkCalibQLook_EXPERT_sig",canvasx,canvasy);
299 c2->SetFillColor(10);
300 tzz->DrawLatex(.01,0.98,fromfile.str().c_str());
301 tzz->DrawLatex(.90,0.98,"SIGMA");
302
303
304 TCanvas *sig=new TCanvas("sig","FTrkCalibQLook_EXPERT_histosig",canvasx,canvasy);
305 sig->SetFillColor(10);
306 tzz->DrawLatex(.01,0.98,fromfile.str().c_str());
307 tzz->DrawLatex(.85,0.97,"Histograms of the sigmas");
308
309
310
311 // draw pads
312 TPad *trkpad1[12],*trkpad2[12],*trkpad3[36]; //pad for histos
313 TPaveText *trkpadtext[12]; //pad for header
314 TH1F *histosig[12]; //histos of sigma
315 TH1F *histoped[12]; //histos of pedestals
316 TH1F *histoasig[12]; //histos of sigma
317 TH1F *histoaped[12]; //histos of pedestals
318
319 TH1F *histosiglad[12][3]; //histos of sigma
320 stringstream title;
321 stringstream hid;
322
323 Float_t posy = 0.95; // up y-coord - top pads
324 Float_t hpad = 0.15; // pad height
325 Float_t posx1=0; // left x-coord - pad column
326 Float_t posx2=0; // right x-coord - pad olumn
327 Float_t posx0=0; // x-coord - column division
328 Float_t wrel = 0.6; // relative x size of first sub-column
329 Float_t marg = 0.004; // margin among pads
330
331
332 for(Int_t n = 0; n<12; n++){
333 if ( (n+1)%2 ) {
334 if(n>1)posy = posy-(marg*2+hpad);
335 posx1 = marg;
336 posx2 = 0.5 - marg;
337 posx0 = 0.5*wrel;
338
339 } else {
340 posx1 = posx1 + 0.5;
341 posx2 = posx2 + 0.5;
342 posx0 = posx0 + 0.5;
343
344 };
345 /* -----------> pad for histograms */
346 trkpad1[n] = new TPad("pad1"," ",posx1,posy-hpad,posx0-marg,posy,18,0,0);
347 trkpad2[n] = new TPad("pad2"," ",posx1,posy-hpad,posx0-marg,posy,18,0,0);
348 /* -----------> pad for header dump */
349 trkpadtext[n] = new TPaveText((posx0+marg),(posy-hpad),posx2,posy);
350 /* -----------> HISTOGRAMS */
351 /* calibration parameters */
352 title<<"DSP "<<n+1;
353 hid<<"h"<<n<<"i"<<i;
354 histosig[n] = new TH1F(hid.str().c_str(),title.str().c_str(),3072,0.5,3072.5);
355 hid.str("");
356 hid<<"hh"<<n<<"i"<<i;
357 histoped[n] = new TH1F(hid.str().c_str(),title.str().c_str(),3072,0.5,3072.5);
358 hid.str("");
359 hid<<"hhh"<<n<<"i"<<i;
360 hid.str("");
361 /* AVERAGE calibration parameters */
362 hid<<"ah"<<n<<"i"<<i;
363 histoasig[n] = new TH1F(hid.str().c_str(),title.str().c_str(),3,0.5,3072.5);
364 hid.str("");
365 hid<<"ahh"<<n<<"i"<<i;
366 histoaped[n] = new TH1F(hid.str().c_str(),title.str().c_str(),3,0.5,3072.5);
367 hid.str("");
368 for(int ii=0;ii<3;ii++){
369 title.str("");
370 title<<"DSP "<<n+1<<" / Lad "<<ii+1;
371 hid<<"hhhh"<<n<<"i"<<i<<"ii"<<ii;
372 histosiglad[n][ii] = new TH1F(hid.str().c_str(),title.str().c_str(),32,-0.5,30.5);
373 hid.str("");
374 }
375 title.str("");
376 }; //end loop on views
377
378 Float_t tposy = 0.95; // up y-coord - top pads
379 Float_t thpad = 0.; // pad height
380 Float_t tposx1=0; // left x-coord - pad column
381 Float_t tposx0=0; // x-coord - column division
382 Float_t twrel = 0.; // relative x size of first sub-column
383 Float_t tmarg = 0.002; // margin among pads
384 thpad = (tposy-tmarg*11)/6;
385 twrel = (1-tmarg*12)/6;
386
387 for(Int_t n = 0; n<36; n++){
388 if ( (n+1)%6==1 ) {
389 if(n>1) tposy = tposy-(tmarg*2+thpad);
390 tposx1 = tmarg;
391 tposx0 = tposx1 + twrel;
392 } else {
393 tposx1 = tposx0 + 2*tmarg;
394 tposx0 = tposx1 + twrel;
395 }
396 trkpad3[n]= new TPad("pad3"," ",tposx1,tposy-thpad,tposx0,tposy,18,0,0);
397 }
398
399
400 stringstream message;
401
402 //--------------------------------
403 //CHECK CALIBRATION procedure
404 //--------------------------------
405 Int_t nn=0;
406 Int_t calok = 0;//BAD
407 for(Int_t n = 0; n<12; n++){
408 if(ctrk.ncalev[n]==0 && ctrk.calfl[n]==0)calok = 1;//GOOD
409
410
411 nn=ctrk.dspnum[n]-1;
412 /*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*
413 *
414 * Write event LEVEL0 report
415 *
416 *.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*/
417 trkpadtext[nn]->SetTextFont(40);
418 trkpadtext[nn]->SetFillColor(33);
419 if(calok==0)trkpadtext[nn]->SetFillColor(2);
420 trkpadtext[nn]->SetTextSize(0.013);
421 trkpadtext[nn]->SetTextAlign(13);
422 trkpadtext[nn]->SetTextColor(1);
423
424 trkpadtext[nn]->AddText(" ");
425 message<<"DSP num --------> "<<ctrk.dspnum[n];
426 trkpadtext[nn]->AddText(message.str().c_str());
427 message.str("");
428 message<<"DAQ mode -------> "<<ctrk.daqmode[n];
429 trkpadtext[nn]->AddText(message.str().c_str());
430 message.str("");
431 message<<"Calibnumber ----> "<<ctrk.calibnum[n];
432 trkpadtext[nn]->AddText(message.str().c_str());
433 message.str("");
434 message<<"ncalib_event ---> "<<ctrk.ncalev[n];
435 trkpadtext[nn]->AddText(message.str().c_str());
436 message.str("");
437 message<<"Cal_flag -------> "<<ctrk.calfl[n];
438 trkpadtext[nn]->AddText(message.str().c_str());
439 message.str("");
440 message<<"crc_hcal -------> "<<ctrk.crc_hc[n];
441 trkpadtext[nn]->AddText(message.str().c_str());
442 message.str("");
443 message<<"crc_cal: l1- ";message.width(5);message<< ctrk.crc_c[n][0];
444 message<<" l2- ";message.width(5);message<< ctrk.crc_c[n][1];
445 message<<" l3- ";message.width(5);message<< ctrk.crc_c[n][2];
446 trkpadtext[nn]->AddText(message.str().c_str());
447 message.str("");
448 message<<"<PED>: l1- "; message.width(5);message<< ctrk.ped1[n];
449 message<<" l2-"; message.width(5);message<< ctrk.ped2[n];
450 message<<" l3-"; message.width(5);message<< ctrk.ped3[n];
451 trkpadtext[nn]->AddText(message.str().c_str());
452 message.str("");
453 message<<"<SIG>: l1- "; message.width(5);message<< ctrk.sig1[n];
454 message<<" l2-"; message.width(5);message<< ctrk.sig2[n];
455 message<<" l3-"; message.width(5);message<< ctrk.sig3[n];
456 trkpadtext[nn]->AddText(message.str().c_str());
457 message.str("");
458 message<<"n.BAD: l1- "; message.width(5);message<<ctrk.nbad1[n];
459 message<<" l2-"; message.width(5);message<< ctrk.nbad2[n];
460 message<<" l3-"; message.width(5);message<< ctrk.nbad3[n];
461 trkpadtext[nn]->AddText(message.str().c_str());
462 message.str("");
463
464 if(calok==0){
465
466 message.str("*** FAILED ***");
467 TText *t=trkpadtext[nn]->AddText(message.str().c_str());
468 t->SetTextSize(0.02);
469 t->SetTextAlign(23);
470 t->SetTextColor(5);
471 };
472
473 /*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*
474 *
475 * Plot histo
476 *
477 *.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*/
478
479
480 /******************************************************/
481 /* fill histos */
482 for(Int_t j = 0; j < 3072; j++){
483 histosig[nn]->Fill((Float_t)j,ctrk.dspsig[nn][j]);
484 histoped[nn]->Fill((Float_t)j,ctrk.dspped[nn][j]);
485 if(j<1024) histosiglad[nn][0]->Fill(ctrk.dspsig[nn][j]);
486 if(j>=1024 && j<2048) histosiglad[nn][1]->Fill(ctrk.dspsig[nn][j]);
487 if(j>=2048 && j<3072) histosiglad[nn][2]->Fill(ctrk.dspsig[nn][j]);
488 };
489 histoasig[nn]->Fill(1,ctrk.sig1[nn]);
490 histoasig[nn]->Fill(1025,ctrk.sig2[nn]);
491 histoasig[nn]->Fill(2049,ctrk.sig3[nn]);
492 histoaped[nn]->Fill(1,ctrk.ped1[nn]);
493 histoaped[nn]->Fill(1025,ctrk.ped2[nn]);
494 histoaped[nn]->Fill(2049,ctrk.ped3[nn]);
495 /******************************************************/
496
497 TLine li;
498 li.SetLineColor(38);
499 li.SetLineStyle(3);
500 li.SetLineWidth(2);
501
502 Float_t maxhist=0;
503 TBox b;
504 /* plot PEDESTAL */
505 c1->cd();
506 trkpadtext[nn]->Draw();
507 trkpad1[nn]->Draw();
508 trkpad1[nn]->cd();
509 trkpad1[nn]->SetFillColor(10);
510 trkpad1[nn]->SetFrameFillColor(10);
511 histoped[nn]->SetStats(kFALSE);
512 histoped[nn]->SetLineColor(1);
513 histoped[nn]->SetFillColor(12);
514 histoped[nn]->SetLineWidth(1);
515 histoped[nn]->GetYaxis()->SetTitle("PED (ADC channels)");
516 histoped[nn]->GetYaxis()->CenterTitle();
517 if((nn+1)%2==1) histoped[nn]->GetYaxis()->SetRangeUser(2200,3200);
518 if((nn+1)%2==0) histoped[nn]->GetYaxis()->SetRangeUser(700,1700);
519 histoaped[nn]->SetLineColor(5);
520 histoaped[nn]->SetLineWidth(1);
521 if(ctrk.good0[0]==1 && ctrk.good0[1]==1){
522 histoped[nn]->Draw("b");
523 if(nn==1){
524 maxhist=histoped[nn]->GetMaximum();
525 b.SetFillColor(6);
526 b.SetFillStyle(3945);
527 b.DrawBox(2944.,700.,3060.,maxhist);
528
529 b.SetFillColor(107);
530 b.SetFillStyle(3954);
531 b.DrawBox(2816.,700.,2944.,maxhist);
532 b.DrawBox(2048.,700.,2176.,maxhist);
533 }
534 else if(nn==6){
535 maxhist=histoped[nn]->GetMaximum();
536 b.SetFillColor(6);
537 b.SetFillStyle(3945);
538 b.DrawBox(2560.,2200.,2816.,maxhist);
539 b.DrawBox(1024.,2200.,1792.,maxhist);
540
541 b.SetFillColor(107);
542 b.SetFillStyle(3954);
543 b.DrawBox(512.,2200.,768.,maxhist);
544 }
545 else if(nn==7){
546 maxhist=histoped[nn]->GetMaximum();
547 b.SetFillColor(107);
548 b.SetFillStyle(3954);
549 b.DrawBox(512.,700.,768.,maxhist);
550 }
551 else if(nn==11){
552 maxhist=histoped[nn]->GetMaximum();
553 b.SetFillColor(6);
554 b.SetFillStyle(3945);
555 b.DrawBox(768.,700.,1024.,maxhist);
556
557 b.SetFillColor(107);
558 b.SetFillStyle(3954);
559 b.DrawBox(0.,700.,512.,maxhist);
560 b.DrawBox(1920.,700.,2048.,maxhist);
561 }
562 }
563 else histoped[nn]->Draw("axis");
564 histoaped[nn]->Draw("same");
565 if((nn+1)%2==1) {
566 li.DrawLine(1024.5,2200,1024.5,3200);
567 li.DrawLine(2048.5,2200,2048.5,3200);
568 }
569 if((nn+1)%2==0) {
570 li.DrawLine(1024.5,700,1024.5,1700);
571 li.DrawLine(2048.5,700,2048.5,1700);
572 }
573
574
575
576 /* plot SIGMA */
577 Float_t max=500.;
578 c2->cd();
579 trkpadtext[nn]->Draw();
580 trkpad2[nn]->SetLogy();
581 trkpad2[nn]->Draw();
582 trkpad2[nn]->cd();
583 trkpad2[nn]->SetFillColor(10);
584 trkpad2[nn]->SetFrameFillColor(10);
585 histosig[nn]->SetStats(kFALSE);
586 histosig[nn]->SetLineColor(1);
587 histosig[nn]->SetFillColor(12);
588 histosig[nn]->SetLineWidth(1);
589 histosig[nn]->SetMaximum(max);
590 histosig[nn]->SetMinimum(0.5);
591 histosig[nn]->GetYaxis()->SetTitle("SIG (ADC channels)");
592 histosig[nn]->GetYaxis()->CenterTitle();
593 histoasig[nn]->SetLineColor(5);
594 histoasig[nn]->SetLineWidth(1);
595 if(ctrk.good0[0]==1 && ctrk.good0[1]==1){
596 histosig[nn]->Draw("b");
597 if(nn==1){
598 maxhist=histosig[nn]->GetMaximum();
599 b.SetFillColor(6);
600 b.SetFillStyle(3945);
601 b.DrawBox(2944.,0.,3060.,maxhist);
602
603 b.SetFillColor(107);
604 b.SetFillStyle(3954);
605 b.DrawBox(2816.,0.,2944.,maxhist);
606 b.DrawBox(2048.,0.,2176.,maxhist);
607 }
608 else if(nn==6){
609 maxhist=histosig[nn]->GetMaximum();
610 b.SetFillColor(6);
611 b.SetFillStyle(3945);
612 b.DrawBox(2560.,0.,2816.,maxhist);
613 b.DrawBox(1024.,0.,1792.,maxhist);
614
615 b.SetFillColor(107);
616 b.SetFillStyle(3954);
617 b.DrawBox(512.,0.,768.,maxhist);
618 }
619 else if(nn==7){
620 maxhist=histosig[nn]->GetMaximum();
621 b.SetFillColor(107);
622 b.SetFillStyle(3954);
623 b.DrawBox(512.,0.,768.,maxhist);
624 }
625 else if(nn==11){
626 maxhist=histosig[nn]->GetMaximum();
627 b.SetFillColor(6);
628 b.SetFillStyle(3945);
629 b.DrawBox(768.,0.,1024.,maxhist);
630
631 b.SetFillColor(107);
632 b.SetFillStyle(3954);
633 b.DrawBox(0.,0.,512.,maxhist);
634 b.DrawBox(1920.,0.,2048.,maxhist);
635 }
636 }
637 else histosig[nn]->Draw("axis");
638 histoasig[nn]->Draw("same");
639 li.DrawLine(1024.5,0,1024.5,max);
640 li.DrawLine(2048.5,0,2048.5,max);
641
642 for(int ii=0;ii<3;ii++){
643 sig->cd();
644 trkpad3[nn*3+ii]->Draw();
645 trkpad3[nn*3+ii]->cd();
646 trkpad3[nn*3+ii]->SetFillColor(10);
647 trkpad3[nn*3+ii]->SetFrameFillColor(10);
648 trkpad3[nn*3+ii]->SetLogy();
649 histosiglad[nn][ii]->SetLineColor(1);
650 histosiglad[nn][ii]->SetFillColor(1);
651 histosiglad[nn][ii]->SetLineWidth(1);
652 histosiglad[nn][ii]->GetXaxis()->SetTitle("SIG (ADC channels)");
653 histosiglad[nn][ii]->GetXaxis()->CenterTitle();
654 histosiglad[nn][ii]->Draw("");
655 }
656
657 };//end loop on views
658 c1->Update();//draw pads in canvas
659 c2->Update();//draw pads in canvas
660 sig->Update();//draw pads in canvas
661
662 stringstream nom1,nom2,nom3;
663
664 if(!strcmp(outfile.Data(),"ps")||!strcmp(outfile.Data(),"pdf")){
665 nom1.str("");
666 nom2.str("");
667 nom3.str("");
668 nom1<<out<<ffile<<"_FTrkCalibQLook_EXPERT-pkt"<<i+1<<".ps(";
669 nom2<<out<<ffile<<"_FTrkCalibQLook_EXPERT-pkt"<<i+1<<".ps";
670 nom3<<out<<ffile<<"_FTrkCalibQLook_EXPERT-pkt"<<i+1<<".ps)";
671 c1->Print(nom1.str().c_str(),"Landscape");
672 c2->Print(nom2.str().c_str(),"Landscape");
673 sig->Print(nom3.str().c_str(),"Landscape");
674
675 if(!strcmp(outfile.Data(),"pdf")){
676 stringstream com;
677 com<<"ps2pdf13 "<<out<<ffile<<"_FTrkCalibQLook_EXPERT-pkt"<<i+1<<".ps "<<out<<ffile<<"_FTrkCalibQlook_EXPERT-pkt"<<i+1<<".pdf";
678 system(com.str().c_str());
679 printf("\n---> ps file converted in pdf format!\n");
680 com.str("");
681 com<<"rm -f "<<out<<ffile<<"_FTrkCalibQLook_EXPERT-pkt"<<i+1<<".ps";
682 system(com.str().c_str());
683 printf("---> ps file removed!\n\n");
684 com.str("");
685 }
686 }
687 else{
688 nom1.str("");
689 nom2.str("");
690 nom3.str("");
691 nom1<<out<<ffile<<"_FTrkCalibQLook_EXPERT-ped-pkt"<<i+1<<"."<<outfile.Data();
692 nom2<<out<<ffile<<"_FTrkCalibQLook_EXPERT-sig-pkt"<<i+1<<"."<<outfile.Data();
693 nom3<<out<<ffile<<"_FTrkCalibQLook_EXPERT-histosig-pkt"<<i+1<<"."<<outfile.Data();
694 c1->Print(nom1.str().c_str());
695 c2->Print(nom2.str().c_str());
696 sig->Print(nom3.str().c_str());
697 }
698
699 if(step==1 && i!=maxevent-1 ){
700 printf("\n Press 1<enter> to continue, 2<enter> to quit.\n");
701 cin>>risposta;
702 if ( i != maxevent-1 ) {
703 if ( risposta == 2 ) {
704 gROOT->Reset();
705 printf("Exiting...\n");
706 return;
707 }
708 }
709 }
710 };//end loop on events
711
712 gROOT->Reset();
713 printf("... end of packets. \n");
714 return;
715 }
716
717

  ViewVC Help
Powered by ViewVC 1.1.23