1 |
/** |
2 |
* FTrkCalibQLook_BASIC.cxx |
3 |
* |
4 |
* autor: D.Fedele |
5 |
* version v1r04 |
6 |
* Parameters: |
7 |
* file - the data file to analyze |
8 |
* fromevent - first event to analyze |
9 |
* toevent - last event to analyze |
10 |
* outdir - total path of output file |
11 |
* outfile - extension of output file (pdf,ps,gif,jpg) |
12 |
* |
13 |
*/ |
14 |
// |
15 |
#include <TLatex.h> |
16 |
#include <TCanvas.h> |
17 |
#include <TGraph.h> |
18 |
#include <TFile.h> |
19 |
#include <TTree.h> |
20 |
#include <TStyle.h> |
21 |
#include <TString.h> |
22 |
// |
23 |
#include <PscuHeader.h> |
24 |
#include <EventHeader.h> |
25 |
#include <CalibTrk1Event.h> |
26 |
#include <CalibTrk2Event.h> |
27 |
// |
28 |
|
29 |
typedef struct caltrk_def{ |
30 |
Int_t good0[2]; |
31 |
Int_t daqmode[12]; |
32 |
Int_t dspnum[12]; |
33 |
Int_t calibnum[12]; |
34 |
Int_t ncalev[12]; |
35 |
Int_t calfl[12]; |
36 |
Int_t ped1[12]; |
37 |
Int_t ped2[12]; |
38 |
Int_t ped3[12]; |
39 |
Int_t sig1[12]; |
40 |
Int_t sig2[12]; |
41 |
Int_t sig3[12]; |
42 |
Int_t nbad1[12]; |
43 |
Int_t nbad2[12]; |
44 |
Int_t nbad3[12]; |
45 |
Float_t dspped[12][3702]; |
46 |
Float_t dspsig[12][3702]; |
47 |
Float_t dspbad[12][3702]; |
48 |
Int_t crc_c[12][3]; |
49 |
Int_t crc_hc[12]; |
50 |
}; |
51 |
|
52 |
|
53 |
void FTrkCalibQLook_BASIC(TString file,Int_t fromevent,Int_t toevent,TString outdir,TString outfile) |
54 |
{ |
55 |
// |
56 |
// obtain information about the data file and select the output dir |
57 |
const string filepath=file.Data(); |
58 |
Int_t dwpos = file.Last('/'); |
59 |
Int_t dwpos1 = file.Last('.'); |
60 |
TString base,ffile ; |
61 |
ffile=file(dwpos+1,dwpos1-(dwpos+1)); |
62 |
if(dwpos>0) base=file(0,dwpos); |
63 |
|
64 |
TString out; |
65 |
if(outdir.Length()==0){ |
66 |
out = base; |
67 |
}else{ |
68 |
out = outdir; |
69 |
}; |
70 |
if(out.Last('/')+1<out.Length()) out+="/"; |
71 |
|
72 |
// |
73 |
// inizialise the variables and open the file |
74 |
struct caltrk_def ctrk; |
75 |
Int_t nevents=0; |
76 |
Int_t minevent = 0; |
77 |
Int_t maxevent = 0; |
78 |
ULong64_t OBT[2]; |
79 |
|
80 |
OBT[0]=0; |
81 |
OBT[1]=0; |
82 |
ctrk.good0[0]=0; |
83 |
ctrk.good0[1]=0; |
84 |
for(Int_t i=0;i<12;i++){ |
85 |
ctrk.daqmode[i]=0; |
86 |
ctrk.dspnum[i]=0; |
87 |
ctrk.calibnum[i]=0; |
88 |
ctrk.ncalev[i]=0; |
89 |
ctrk.calfl[i]=0; |
90 |
ctrk.ped1[i]=0; |
91 |
ctrk.ped2[i]=0; |
92 |
ctrk.ped3[i]=0; |
93 |
ctrk.sig1[i]=0; |
94 |
ctrk.sig2[i]=0; |
95 |
ctrk.sig3[i]=0; |
96 |
ctrk.nbad1[i]=0; |
97 |
ctrk.nbad2[i]=0; |
98 |
ctrk.nbad3[i]=0; |
99 |
ctrk.crc_hc[i]=0; |
100 |
ctrk.crc_c[i][0]=0; |
101 |
ctrk.crc_c[i][1]=0; |
102 |
ctrk.crc_c[i][2]=0; |
103 |
for(Int_t iii=0;iii<3072;iii++){ |
104 |
ctrk.dspped[i][iii]=0; |
105 |
ctrk.dspsig[i][iii]=0; |
106 |
ctrk.dspbad[i][iii]=0; |
107 |
} |
108 |
} |
109 |
|
110 |
pamela::EventHeader *eh1=0,*eh2=0; |
111 |
pamela::PscuHeader *ph1=0,*ph2=0; |
112 |
pamela::CalibTrk1Event *trk1 = 0; |
113 |
pamela::CalibTrk2Event *trk2 = 0; |
114 |
|
115 |
TFile *datafile = new TFile(file); |
116 |
if ( !datafile ){ |
117 |
printf("No data file, exiting...\n"); |
118 |
return; |
119 |
}; |
120 |
TTree *otr1,*otr2; |
121 |
|
122 |
otr1 = (TTree*)datafile->Get("CalibTrk1"); |
123 |
otr1->SetBranchAddress("CalibTrk1", &trk1); |
124 |
otr1->SetBranchAddress("Header",&eh1); |
125 |
otr2 = (TTree*)datafile->Get("CalibTrk2"); |
126 |
otr2->SetBranchAddress("CalibTrk2", &trk2); |
127 |
otr2->SetBranchAddress("Header",&eh2); |
128 |
|
129 |
|
130 |
if(otr1->GetEntries()==otr2->GetEntries()) |
131 |
nevents = otr1->GetEntries(); |
132 |
else{ |
133 |
printf("WARNING: CalibTrk1 entries is different from CalibTrk2 entries"); |
134 |
return;} |
135 |
|
136 |
if (nevents<=0) { |
137 |
datafile->Close(); |
138 |
printf("No calibration packets found, exiting...\n"); |
139 |
return; |
140 |
}; |
141 |
printf("Number of calibration packets: %i\n",nevents); |
142 |
|
143 |
if ( fromevent > toevent && toevent ){ |
144 |
printf("It must be fromevent < toevent \n"); |
145 |
return; |
146 |
}; |
147 |
|
148 |
if ( fromevent > nevents || fromevent < 0 ) { |
149 |
printf("You can choose fromevent between 0 (all) and %i \n",nevents); |
150 |
return; |
151 |
}; |
152 |
|
153 |
if ( toevent > nevents || toevent < 0 ) { |
154 |
printf("You can choose toevent between 0 (all) and %i \n",nevents); |
155 |
return; |
156 |
}; |
157 |
if ( fromevent == 0 ) { |
158 |
minevent = 0; |
159 |
maxevent = nevents; |
160 |
} else { |
161 |
minevent = fromevent - 1; |
162 |
if ( toevent > 0 ){ |
163 |
maxevent = toevent; |
164 |
} else if (toevent > nevents) { |
165 |
maxevent = nevents; |
166 |
} else { |
167 |
maxevent = fromevent; |
168 |
}; |
169 |
}; |
170 |
|
171 |
// |
172 |
// other variables definitions |
173 |
stringstream fromfile,rep,tit; |
174 |
fromfile<<"FTrkCalibQLook_BASIC File: "<<ffile; |
175 |
|
176 |
gStyle->SetLabelSize(0.08,"x"); |
177 |
gStyle->SetLabelSize(0.08,"y"); |
178 |
gStyle->SetFillColor(0); |
179 |
gStyle->SetTitleFillColor(0); |
180 |
gStyle->SetTitleFontSize(0.1); |
181 |
gStyle->SetTitleOffset(0.8,"y"); |
182 |
gStyle->SetTitleOffset(1.,"x"); |
183 |
gStyle->SetTitleSize(0.06,"y"); |
184 |
gStyle->SetTitleSize(0.06,"x"); |
185 |
gStyle->SetOptStat(0); |
186 |
|
187 |
TLatex *tzz=new TLatex(); |
188 |
tzz->SetTextFont(32); |
189 |
tzz->SetTextColor(1); |
190 |
tzz->SetTextAlign(12); |
191 |
tzz->SetTextSize(0.02); |
192 |
|
193 |
Int_t canvasx=900; |
194 |
Int_t canvasy=1200; |
195 |
|
196 |
Int_t ndsp =0,alarm=0; |
197 |
Float_t pedav[12][12],pedavtemp[12][12],sigav[12][12],sigavtemp[12][12]; |
198 |
Int_t flpedav[12][12],flsigav[12][12]; |
199 |
Float_t siglimsup[12][12],sigliminf[12][12],pedlimsup[12][12],pedliminf[12][12]; |
200 |
|
201 |
// |
202 |
// inizialize the limits for simga and pedestall |
203 |
for(Int_t i=0;i<12;i++){ |
204 |
for(Int_t ii=0;ii<12;ii++){ |
205 |
siglimsup[i][ii]=30.; |
206 |
sigliminf[i][ii]=1.5; |
207 |
if(!(i%2)){ |
208 |
pedlimsup[i][ii]=3700.; |
209 |
pedliminf[i][ii]=1700.; |
210 |
} |
211 |
else{ |
212 |
pedlimsup[i][ii]=2200.; |
213 |
pedliminf[i][ii]=200.; |
214 |
} |
215 |
} |
216 |
} |
217 |
|
218 |
// |
219 |
// count of possible alarm to set the number of output pages |
220 |
alarm+=(maxevent-minevent)/9 +1; |
221 |
for (Int_t i = minevent; i < maxevent; i++){ |
222 |
otr1->GetEntry(i); |
223 |
otr2->GetEntry(i); |
224 |
for(Int_t m=0;m<6;m++){ |
225 |
for(Int_t mm=0;mm<3;mm++){ |
226 |
if(trk1->crc_cal[m][mm]!=0) alarm+=1; |
227 |
if(trk2->crc_cal[m][mm]!=0 ) alarm+=1; |
228 |
} |
229 |
if(trk1->crc_hcal[m]!=0 ) alarm+=1; |
230 |
if(trk1->cal_flag[m]!=0 ) alarm+=1; |
231 |
if(trk2->crc_hcal[m]!=0) alarm+=1; |
232 |
if(trk2->cal_flag[m]!=0) alarm+=1; |
233 |
} |
234 |
} |
235 |
const Int_t cnum=alarm/30 + 100; |
236 |
Int_t flcanvas=1; |
237 |
Float_t spacep=1.5,space[cnum]; |
238 |
TCanvas *c[cnum]; |
239 |
// |
240 |
// create output canvas |
241 |
for(Int_t i=0;i<cnum;i++){ |
242 |
space[i]=96.0; |
243 |
rep.str(""); |
244 |
tit.str(""); |
245 |
tit<<"c"<<i; |
246 |
rep<<"FTrkCalibQLook_BASIC_pag"<<i+1; |
247 |
c[i]=new TCanvas(tit.str().c_str(),rep.str().c_str(),canvasx,canvasy); |
248 |
c[i]->Range(0,0,100,100); |
249 |
c[i]->SetFillColor(10); |
250 |
tzz->DrawLatex(1,98.5,fromfile.str().c_str()); |
251 |
rep.str(""); |
252 |
rep<<"CALIBRATION REPORT pag"<<i+1; |
253 |
tzz->DrawLatex(70,98.5,rep.str().c_str()); |
254 |
rep.str(""); |
255 |
} |
256 |
|
257 |
|
258 |
//********************************************************************** |
259 |
// |
260 |
// LOOP OVER EVENTS |
261 |
// |
262 |
//********************************************************************** |
263 |
|
264 |
Int_t wc=0; |
265 |
printf("\n Scan of calibration packets from %i to %i ... \n",minevent+1,maxevent); |
266 |
for (Int_t i = minevent; i < maxevent; i++){ |
267 |
|
268 |
otr1->GetEntry(i); |
269 |
otr2->GetEntry(i); |
270 |
|
271 |
ctrk.good0[0]=trk1->good0; |
272 |
ctrk.good0[1]=trk2->good0; |
273 |
for (Int_t m = 0; m < 6; m++){ |
274 |
ph1 = eh1->GetPscuHeader(); |
275 |
OBT[0]= ph1->GetOrbitalTime(); |
276 |
ctrk.daqmode[trk1->DSPnumber[m]-1]=trk1->DAQmode[m]; |
277 |
ctrk.dspnum[trk1->DSPnumber[m]-1]=trk1->DSPnumber[m]; |
278 |
ctrk.calibnum[trk1->DSPnumber[m]-1]=trk1->calibnumber[m]; |
279 |
ctrk.ncalev[trk1->DSPnumber[m]-1]=trk1->ncalib_event[m]; |
280 |
ctrk.ped1[trk1->DSPnumber[m]-1]=trk1->ped_l1[m]; |
281 |
ctrk.ped2[trk1->DSPnumber[m]-1]=trk1->ped_l2[m]; |
282 |
ctrk.ped3[trk1->DSPnumber[m]-1]=trk1->ped_l3[m]; |
283 |
ctrk.sig1[trk1->DSPnumber[m]-1]=trk1->sig_l1[m]; |
284 |
ctrk.sig2[trk1->DSPnumber[m]-1]=trk1->sig_l2[m]; |
285 |
ctrk.sig3[trk1->DSPnumber[m]-1]=trk1->sig_l3[m]; |
286 |
ctrk.nbad1[trk1->DSPnumber[m]-1]=trk1->nbad_l1[m]; |
287 |
ctrk.nbad2[trk1->DSPnumber[m]-1]=trk1->nbad_l2[m]; |
288 |
ctrk.nbad3[trk1->DSPnumber[m]-1]=trk1->nbad_l3[m]; |
289 |
ctrk.calfl[trk1->DSPnumber[m]-1]=trk1->cal_flag[m]; |
290 |
ctrk.crc_c[trk1->DSPnumber[m]-1][0]=trk1->crc_cal[m][0]; |
291 |
ctrk.crc_c[trk1->DSPnumber[m]-1][1]=trk1->crc_cal[m][1]; |
292 |
ctrk.crc_c[trk1->DSPnumber[m]-1][2]=trk1->crc_cal[m][2]; |
293 |
ctrk.crc_hc[trk1->DSPnumber[m]-1]=trk1->crc_hcal[m]; |
294 |
for (Int_t j = 0; j < 3072; j++){ |
295 |
ctrk.dspped[trk1->DSPnumber[m]-1][j]=trk1->DSPped_par[m][j]; |
296 |
ctrk.dspsig[trk1->DSPnumber[m]-1][j]=trk1->DSPsig_par[m][j]; |
297 |
ctrk.dspbad[trk1->DSPnumber[m]-1][j]=trk1->DSPbad_par[m][j]; |
298 |
} |
299 |
ph2 = eh2->GetPscuHeader(); |
300 |
OBT[1]= ph2->GetOrbitalTime(); |
301 |
ctrk.daqmode[trk2->DSPnumber[m]-1]=trk2->DAQmode[m]; |
302 |
ctrk.dspnum[trk2->DSPnumber[m]-1]=trk2->DSPnumber[m]; |
303 |
ctrk.calibnum[trk2->DSPnumber[m]-1]=trk2->calibnumber[m]; |
304 |
ctrk.ncalev[trk2->DSPnumber[m]-1]=trk2->ncalib_event[m]; |
305 |
ctrk.ped1[trk2->DSPnumber[m]-1]=trk2->ped_l1[m]; |
306 |
ctrk.ped2[trk2->DSPnumber[m]-1]=trk2->ped_l2[m]; |
307 |
ctrk.ped3[trk2->DSPnumber[m]-1]=trk2->ped_l3[m]; |
308 |
ctrk.sig1[trk2->DSPnumber[m]-1]=trk2->sig_l1[m]; |
309 |
ctrk.sig2[trk2->DSPnumber[m]-1]=trk2->sig_l2[m]; |
310 |
ctrk.sig3[trk2->DSPnumber[m]-1]=trk2->sig_l3[m]; |
311 |
ctrk.nbad1[trk2->DSPnumber[m]-1]=trk2->nbad_l1[m]; |
312 |
ctrk.nbad2[trk2->DSPnumber[m]-1]=trk2->nbad_l2[m]; |
313 |
ctrk.nbad3[trk2->DSPnumber[m]-1]=trk2->nbad_l3[m]; |
314 |
ctrk.calfl[trk2->DSPnumber[m]-1]=trk2->cal_flag[m]; |
315 |
ctrk.crc_c[trk1->DSPnumber[m]-1][0]=trk2->crc_cal[m][0]; |
316 |
ctrk.crc_c[trk1->DSPnumber[m]-1][1]=trk2->crc_cal[m][1]; |
317 |
ctrk.crc_c[trk1->DSPnumber[m]-1][2]=trk2->crc_cal[m][2]; |
318 |
ctrk.crc_hc[trk1->DSPnumber[m]-1]=trk2->crc_hcal[m]; |
319 |
for (Int_t j = 0; j < 3072; j++){ |
320 |
ctrk.dspped[trk2->DSPnumber[m]-1][j]=trk2->DSPped_par[m][j]; |
321 |
ctrk.dspsig[trk2->DSPnumber[m]-1][j]=trk2->DSPsig_par[m][j]; |
322 |
ctrk.dspbad[trk2->DSPnumber[m]-1][j]=trk2->DSPbad_par[m][j]; |
323 |
} |
324 |
} |
325 |
|
326 |
for(Int_t n = 0; n<12; n++){ |
327 |
for(Int_t nm = 0; nm<12; nm++){ |
328 |
pedav[n][nm]=0; |
329 |
pedavtemp[n][nm]=0; |
330 |
sigav[n][nm]=0; |
331 |
sigavtemp[n][nm]=0; |
332 |
flpedav[n][nm]=0; |
333 |
flsigav[n][nm]=0; |
334 |
} |
335 |
} |
336 |
|
337 |
Int_t nn,ok=0; |
338 |
|
339 |
// |
340 |
// write warning if it occur |
341 |
for(Int_t n = 0; n<12; n++){ |
342 |
|
343 |
ndsp = ctrk.dspnum[n]; |
344 |
nn = ndsp-1; |
345 |
|
346 |
for(Int_t iii=0;iii<3;iii++){ |
347 |
if(ctrk.crc_c[nn][iii]!=0){ |
348 |
ok=1; |
349 |
if(space[wc]<=3){ |
350 |
wc+=1; |
351 |
flcanvas+=1; |
352 |
} |
353 |
c[wc]->cd(); |
354 |
tzz->SetTextFont(40); |
355 |
tzz->SetTextSize(0.02); |
356 |
tzz->SetTextAlign(13); |
357 |
tzz->SetTextColor(2); |
358 |
rep<<"***************************************************************************************************************************"; |
359 |
tzz->DrawLatex(2,space[wc],rep.str().c_str()); |
360 |
rep.str(""); |
361 |
space[wc]-=spacep; |
362 |
rep<<" ERROR >>> CALIBRATION pkt "<<i+1<<" -->CalibTrk"<<(nn+1)%2+1<<" at OBT: "<<OBT[(nn+1)%2]<<" --> crc_cal["<<nn+1<<"]["<<iii+1<<"]= "<<ctrk.crc_c[nn][iii]; |
363 |
tzz->DrawLatex(2,space[wc],rep.str().c_str()); |
364 |
rep.str(""); |
365 |
space[wc]-=spacep; |
366 |
} |
367 |
} |
368 |
|
369 |
if(ctrk.crc_hc[nn]!=0){ |
370 |
ok=1; |
371 |
if(space[wc]<=3){ |
372 |
wc+=1; |
373 |
flcanvas+=1; |
374 |
} |
375 |
c[wc]->cd(); |
376 |
tzz->SetTextFont(40); |
377 |
tzz->SetTextSize(0.02); |
378 |
tzz->SetTextAlign(13); |
379 |
tzz->SetTextColor(2); |
380 |
rep<<"***************************************************************************************************************************"; |
381 |
tzz->DrawLatex(2,space[wc],rep.str().c_str()); |
382 |
rep.str(""); |
383 |
space[wc]-=spacep; |
384 |
rep<<" ERROR >>> CALIBRATION pkt "<<i+1<<" -->CalibTrk"<<(nn+1)%2+1<<" at OBT: "<<OBT[(nn+1)%2]<<" --> crc_hcal["<<nn+1<<"]= "<<ctrk.crc_hc[nn]; |
385 |
tzz->DrawLatex(2,space[wc],rep.str().c_str()); |
386 |
rep.str(""); |
387 |
space[wc]-=spacep; |
388 |
} |
389 |
|
390 |
if(ctrk.calfl[nn]!=0){ |
391 |
ok=1; |
392 |
if(space[wc]<=3){ |
393 |
wc+=1; |
394 |
flcanvas+=1; |
395 |
} |
396 |
c[wc]->cd(); |
397 |
tzz->SetTextFont(40); |
398 |
tzz->SetTextSize(0.02); |
399 |
tzz->SetTextAlign(13); |
400 |
tzz->SetTextColor(2); |
401 |
rep<<"***************************************************************************************************************************"; |
402 |
tzz->DrawLatex(2,space[wc],rep.str().c_str()); |
403 |
rep.str(""); |
404 |
space[wc]-=spacep; |
405 |
rep<<" ERROR >>> CALIBRATION pkt "<<i+1<<" -->CalibTrk"<<(nn+1)%2+1<<" at OBT: "<<OBT[(nn+1)%2]<<" --> cal_flag["<<nn+1<<"]= "<<ctrk.calfl[nn]; |
406 |
tzz->DrawLatex(2,space[wc],rep.str().c_str()); |
407 |
rep.str(""); |
408 |
space[wc]-=spacep; |
409 |
} |
410 |
|
411 |
if(ctrk.good0[0]==1 && ctrk.good0[1]==1){ |
412 |
// |
413 |
// evaluate the mean value of the sigma and pedestal |
414 |
for(Int_t j = 0; j < 3072; j++){ |
415 |
sigavtemp[nn][j/256]+=ctrk.dspsig[nn][j]; |
416 |
pedavtemp[nn][j/256]+=ctrk.dspped[nn][j]; |
417 |
} |
418 |
|
419 |
for(Int_t ii=0;ii<12;ii++){ |
420 |
pedav[nn][ii]=pedavtemp[nn][ii]/256; |
421 |
sigav[nn][ii]=sigavtemp[nn][ii]/256; |
422 |
|
423 |
if(pedav[nn][ii]>pedlimsup[nn][ii] || pedav[nn][ii]<pedliminf[nn][ii]) flpedav[nn][ii]=1; |
424 |
if(sigav[nn][ii]>siglimsup[nn][ii] || sigav[nn][ii]<sigliminf[nn][ii]) flsigav[nn][ii]=1; |
425 |
|
426 |
if((nn==1 && ii==11)||(nn==6 && ii==2)||(nn==6 && ii==4)||(nn==6 && ii==5)||(nn==6 && ii==6)||(nn==6 && ii==10)||(nn==11 && ii==3)) |
427 |
continue; |
428 |
else{ |
429 |
if(flpedav[nn][ii]==1){ |
430 |
ok=1; |
431 |
if(space[wc]<=3){ |
432 |
wc+=1; |
433 |
flcanvas+=1; |
434 |
} |
435 |
c[wc]->cd(); |
436 |
tzz->SetTextFont(40); |
437 |
tzz->SetTextSize(0.02); |
438 |
tzz->SetTextAlign(13); |
439 |
tzz->SetTextColor(50); |
440 |
rep<<"********************************************************************************************************************************"; |
441 |
tzz->DrawLatex(2,space[wc],rep.str().c_str()); |
442 |
rep.str(""); |
443 |
space[wc]-=spacep; |
444 |
rep<<"WARNING >>> CALIBRATION pkt "<<i+1<<" -->CalibTrk"<<(nn+1)%2+1<<" at OBT: "<<OBT[(nn+1)%2]<<"-->DSP "<<nn+1<<" -VA1 "<<2*ii+1<<"-"<<2*ii+2<<" --> <PED>= "<<pedav[nn][ii]; |
445 |
tzz->DrawLatex(2,space[wc],rep.str().c_str()); |
446 |
rep.str(""); |
447 |
space[wc]-=spacep; |
448 |
} |
449 |
|
450 |
if(flsigav[nn][ii]==1){ |
451 |
ok=1; |
452 |
if(space[wc]<=3){ |
453 |
wc+=1; |
454 |
flcanvas+=1; |
455 |
} |
456 |
c[wc]->cd(); |
457 |
tzz->SetTextFont(40); |
458 |
tzz->SetTextSize(0.02); |
459 |
tzz->SetTextAlign(13); |
460 |
tzz->SetTextColor(50); |
461 |
rep<<"********************************************************************************************************************************"; |
462 |
tzz->DrawLatex(2,space[wc],rep.str().c_str()); |
463 |
rep.str(""); |
464 |
space[wc]-=spacep; |
465 |
rep<<"WARNING >>> CALIBRATION pkt "<<i+1<<" -->CalibTrk"<<(nn+1)%2+1<<" at OBT: "<<OBT[(nn+1)%2]<<"-->DSP "<<nn+1<<" -VA1 "<<2*ii+1<<"-"<<2*ii+2<<" --> <SIG>= "<<sigav[nn][ii]; |
466 |
tzz->DrawLatex(2,space[wc],rep.str().c_str()); |
467 |
rep.str(""); |
468 |
space[wc]-=spacep; |
469 |
} |
470 |
} |
471 |
} |
472 |
} |
473 |
} |
474 |
if(ok==0 && ctrk.good0[0]==1 && ctrk.good0[1]==1){ |
475 |
if(space[wc]<=10){ |
476 |
wc+=1; |
477 |
flcanvas+=1; |
478 |
} |
479 |
c[wc]->cd(); |
480 |
tzz->SetTextFont(40); |
481 |
tzz->SetTextSize(0.03); |
482 |
tzz->SetTextAlign(13); |
483 |
tzz->SetTextColor(1); |
484 |
rep<<"*********************************** CALIBRATION pkt "<<i+1<<"********************************"; |
485 |
tzz->DrawLatex(2,space[wc],rep.str().c_str()); |
486 |
rep.str(""); |
487 |
space[wc]-=spacep+1; |
488 |
rep<<">>>>>>> CalibTrk1 at OBT: "<<OBT[0]<<" ---------> OK "; |
489 |
tzz->DrawLatex(2,space[wc],rep.str().c_str()); |
490 |
rep.str(""); |
491 |
space[wc]-=spacep+1; |
492 |
rep<<">>>>>>> CalibTrk2 at OBT: "<<OBT[1]<<" ---------> OK "; |
493 |
tzz->DrawLatex(2,space[wc],rep.str().c_str()); |
494 |
rep.str(""); |
495 |
space[wc]-=spacep+1; |
496 |
rep<<"*******************************************************************************************"; |
497 |
tzz->DrawLatex(2,space[wc],rep.str().c_str()); |
498 |
rep.str(""); |
499 |
space[wc]-=spacep+1; |
500 |
} |
501 |
};//end loop on events |
502 |
|
503 |
printf("... end of packets. \n"); |
504 |
// |
505 |
// Save output Files |
506 |
stringstream nom1,nom2,nom3; |
507 |
|
508 |
for(Int_t i=0;i<cnum;i++) |
509 |
c[i]->Update();//draw pads in canvas |
510 |
|
511 |
for(Int_t fl=0;fl<flcanvas;fl++){ |
512 |
if(flcanvas==1){ |
513 |
nom1.str(""); |
514 |
nom1<<ffile<<"_FTrkCalibQLook_BASIC."<<outfile.Data(); |
515 |
c[fl]->Print(out+nom1.str().c_str()); |
516 |
nom1.str(""); |
517 |
} |
518 |
|
519 |
if(flcanvas>=2){ |
520 |
if(!strcmp(outfile.Data(),"ps")||!strcmp(outfile.Data(),"pdf")){ |
521 |
nom1.str(""); |
522 |
nom2.str(""); |
523 |
nom3.str(""); |
524 |
nom1<<ffile<<"_FTrkCalibQLook_BASIC.ps("; |
525 |
nom2<<ffile<<"_FTrkCalibQLook_BASIC.ps"; |
526 |
nom3<<ffile<<"_FTrkCalibQLook_BASIC.ps)"; |
527 |
if(fl==0) c[fl]->Print(out+nom1.str().c_str(),"portrait"); |
528 |
else if(fl==flcanvas-1) c[fl]->Print(out+nom3.str().c_str(),"portrait"); |
529 |
else c[fl]->Print(out+nom2.str().c_str(),"portrait"); |
530 |
|
531 |
} |
532 |
else{ |
533 |
nom1.str(""); |
534 |
nom1<<ffile<<"_FTrkCalibQLook_BASIC-pag"<<fl+1<<"."<<outfile.Data(); |
535 |
c[fl]->Print(out+nom1.str().c_str()); |
536 |
} |
537 |
} |
538 |
} |
539 |
|
540 |
if(!strcmp(outfile.Data(),"pdf")&&flcanvas>=2){ |
541 |
stringstream com; |
542 |
com<<"ps2pdf13 "<<out<<ffile<<"_FTrkCalibQLook_BASIC.ps "<<out<<ffile<<"_FTrkCalibQLook_BASIC.pdf"; |
543 |
system(com.str().c_str()); |
544 |
printf("\n---> ps file converted in pdf format!\n"); |
545 |
com.str(""); |
546 |
com<<"rm -f "<<out<<ffile<<"_FTrkCalibQLook_BASIC.ps"; |
547 |
system(com.str().c_str()); |
548 |
printf("---> ps file removed!\n\n"); |
549 |
com.str(""); |
550 |
} |
551 |
|
552 |
|
553 |
gROOT->Reset(); |
554 |
return; |
555 |
|
556 |
} |