1 |
/** |
2 |
* FTrkScanQlook_EXPERT.cxx |
3 |
* |
4 |
* autor: D.Fedele |
5 |
* version v2r00 |
6 |
* Parameters: |
7 |
* file - the path to the root file to analyze |
8 |
* outdir - total path of output file |
9 |
* event - the number of the single event to analyze |
10 |
* va1 - the number of the single va1 to analyze (dsp*100+va1) |
11 |
* outfile - extension of output file (pdf,ps,gif,jpg) |
12 |
* |
13 |
*/ |
14 |
// |
15 |
#include <iostream> |
16 |
#include <sstream> |
17 |
// |
18 |
#include <TPaveText.h> |
19 |
#include <TLatex.h> |
20 |
#include <TCanvas.h> |
21 |
#include <TGraph.h> |
22 |
#include <TFile.h> |
23 |
#include <TTree.h> |
24 |
#include <TStyle.h> |
25 |
#include <TString.h> |
26 |
// |
27 |
#include <physics/tracker/TrackerEvent.h> |
28 |
#include <PscuHeader.h> |
29 |
#include <EventHeader.h> |
30 |
#include <RunHeaderEvent.h> |
31 |
// |
32 |
|
33 |
using namespace std; |
34 |
|
35 |
typedef struct trkword{ |
36 |
int type; |
37 |
int decode; |
38 |
}; |
39 |
|
40 |
/*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-* |
41 |
* |
42 |
| The function "datadecode" decodes the tracker words. |
43 |
* |
44 |
| Tracker words are of three types: |
45 |
* - ADC data |
46 |
| - strip address |
47 |
* - end of ladder |
48 |
| |
49 |
* The function returns a struct variable (trkword) |
50 |
| that contains two quantities, type and decode, which are assigned |
51 |
* the following values: |
52 |
| |
53 |
* type decode |
54 |
| ---------------------------------------------------- |
55 |
* -1 error 0 |
56 |
| 0 data 0-4095 ADC values |
57 |
* 1 address 0-1023 strip address (relative to ladder) |
58 |
| 2 end-of-ladder 1,2,3 ladder number (compressed acq mode) |
59 |
* 4,5,6 ladder number + 3 (full acq mode) |
60 |
| |
61 |
*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*/ |
62 |
trkword datadecode(int word){ |
63 |
int type =0; |
64 |
int data =0; |
65 |
int nodata = word&0xf000; |
66 |
int zero = 0; |
67 |
|
68 |
trkword thisword; |
69 |
|
70 |
switch (nodata>>12){ |
71 |
|
72 |
case 0: |
73 |
|
74 |
thisword.decode = word; |
75 |
thisword.type = 0; |
76 |
// cout << thisword.decode << "\n"; |
77 |
return (thisword); //>>>>> ADC data (0) |
78 |
|
79 |
case 1: |
80 |
|
81 |
type = word&0xC00; |
82 |
data = word&0x3ff; |
83 |
|
84 |
switch(type>>10){ |
85 |
|
86 |
case 0: |
87 |
thisword.decode = data; |
88 |
thisword.type = 1; //>>> address (1) |
89 |
return (thisword); |
90 |
|
91 |
case 2: |
92 |
// if(data>=4)data = data-3; |
93 |
if(data>6){ |
94 |
printf("Error on data \n"); |
95 |
thisword.decode = zero; |
96 |
thisword.type = -1; |
97 |
return (thisword); //>>>>> error (-1) |
98 |
} |
99 |
thisword.decode = data; |
100 |
thisword.type = 2; |
101 |
return (thisword); //>>> end-of-ladder |
102 |
|
103 |
default: |
104 |
printf("Error on data \n"); |
105 |
thisword.decode = zero; |
106 |
thisword.type = -1; |
107 |
return (thisword); //>>>>> error (-1) |
108 |
|
109 |
} |
110 |
|
111 |
default: |
112 |
|
113 |
printf("Error on data \n"); |
114 |
thisword.decode = zero; |
115 |
thisword.type = -1; |
116 |
return (thisword); //>>>>> error (-1) |
117 |
} |
118 |
} |
119 |
|
120 |
|
121 |
void FTrkScanQLook_EXPERT(TString file, TString outdir,Int_t event, Int_t va1, TString outfile) |
122 |
{ |
123 |
|
124 |
// |
125 |
// obtain information about the data file and select the output file |
126 |
Int_t dwpos = file.Last('/'); |
127 |
Int_t dwpos1 = file.Last('.'); |
128 |
TString base,ffile ; |
129 |
ffile=file(dwpos+1,dwpos1-(dwpos+1)); |
130 |
if(dwpos>0) base=file(0,dwpos); |
131 |
|
132 |
TString out; |
133 |
if(outdir.Length()==0){ |
134 |
out = base; |
135 |
}else{ |
136 |
out = outdir; |
137 |
} |
138 |
if(out.Last('/')+1<out.Length()) out+="/"; |
139 |
|
140 |
pamela::tracker::TrackerEvent *trk=0; |
141 |
pamela::EventHeader *eh=0,*eH=0; |
142 |
pamela::PscuHeader *ph=0,*pH=0; |
143 |
pamela::RunHeaderEvent *reh=0; |
144 |
pamela::EventCounter *cod=0; |
145 |
|
146 |
pamela::PacketType *pctp=0; |
147 |
// open files |
148 |
TFile *trackerFile = new TFile(file); |
149 |
if ( !trackerFile ){ |
150 |
trackerFile->Close(); |
151 |
printf("No tracker file! \n"); |
152 |
return; |
153 |
} |
154 |
|
155 |
//Takes the tree and branches |
156 |
TTree *tr = (TTree*)trackerFile->Get("Physics"); |
157 |
tr->SetBranchAddress("Tracker",&trk); |
158 |
tr->SetBranchAddress("Header",&eh); |
159 |
|
160 |
TTree *otr = (TTree*)trackerFile->Get("RunHeader"); |
161 |
otr->SetBranchAddress("Header",&eH); |
162 |
otr->SetBranchAddress("RunHeader",&reh); |
163 |
|
164 |
// Define variables |
165 |
Int_t nevents = tr->GetEntries(); |
166 |
Int_t neventH = otr->GetEntries(); |
167 |
if ( nevents <= 0 ) { |
168 |
trackerFile->Close(); |
169 |
printf("The file is empty, exiting...\n"); |
170 |
return; |
171 |
} |
172 |
|
173 |
|
174 |
printf("\n Number of Entries: %d\n",nevents); |
175 |
printf(" Number of Header Entries: %d\n",neventH); |
176 |
|
177 |
Long64_t obt=0; |
178 |
Int_t eve,cin=0;//,p=0; |
179 |
// if(event<0) Int_t eve[abs(event)]; |
180 |
TString cal=""; |
181 |
|
182 |
eve=3; |
183 |
//ev[0]=3; |
184 |
// ev[1]=4; |
185 |
for(Int_t i=0;i<neventH;i++){ |
186 |
otr->GetEntry(i); |
187 |
pH = eH->GetPscuHeader(); |
188 |
if(reh->TRK_CALIB_USED!=104){ |
189 |
obt = pH->GetOrbitalTime(); |
190 |
cal="Event with online calibration"; |
191 |
break; |
192 |
} |
193 |
if(i==neventH-1){ |
194 |
cal="***** ONLINE CALIBRATION NOT FOUND IN THIS FILE *****"; |
195 |
eve=2; |
196 |
// ev[0]=2; |
197 |
//ev[1]=3; |
198 |
} |
199 |
} |
200 |
if(eve==3){ |
201 |
for(Int_t i=0;i<nevents;i++){ |
202 |
tr->GetEntry(i); |
203 |
ph = eh->GetPscuHeader(); |
204 |
cod = eh->GetCounter(); |
205 |
if(i==0) cin=cod->Get(pctp->CalibTrk1); |
206 |
if(reh->TRK_CALIB_USED==104) continue; |
207 |
if(event<0 && cod->Get(pctp->CalibTrk1)==cin+1){ |
208 |
eve=i+3; |
209 |
break; |
210 |
} |
211 |
else if(event>=0 && ph->GetOrbitalTime()>obt){ |
212 |
eve=i+3; |
213 |
//ev[0]=i+3; |
214 |
//ev[1]=i+4; |
215 |
break; |
216 |
} |
217 |
} |
218 |
} |
219 |
|
220 |
int tot=2; |
221 |
if(event<0) tot=abs(event); |
222 |
|
223 |
TH1F *histomax[12][tot]; //histos of max signals |
224 |
TH1F *histocomp[12][tot]; //histos of compressed data |
225 |
TH1F *histofull[12][tot]; //histos of full data |
226 |
TCanvas *c1[tot]; |
227 |
TH1F *histova[tot]; //histos of va1 signals |
228 |
TCanvas *cva[tot]; |
229 |
|
230 |
|
231 |
TLatex *t=new TLatex(); |
232 |
t->SetTextFont(32); |
233 |
t->SetTextColor(1); |
234 |
t->SetTextAlign(12); |
235 |
t->SetTextSize(0.02); |
236 |
|
237 |
|
238 |
for(Int_t e=0;e<tot;e++){ |
239 |
if(event<=0) |
240 |
event=eve; |
241 |
else { |
242 |
event=event+e; |
243 |
if(event>eve-3 && eve>2) |
244 |
cal="Event with online calibration"; |
245 |
else |
246 |
cal="***** ONLINE CALIBRATION NOT FOUND IN THIS FILE *****"; |
247 |
} |
248 |
printf("Scan of Entry %d\n",event); |
249 |
|
250 |
tr->GetEntry(event); |
251 |
//============================================================================ |
252 |
|
253 |
gStyle->SetLabelSize(0.06,"x"); |
254 |
gStyle->SetLabelSize(0.06,"y"); |
255 |
//gStyle->SetTitleFontSize(0.1); |
256 |
gStyle->SetFillColor(10); |
257 |
gStyle->SetTitleFillColor(10); |
258 |
gStyle->SetTitleOffset(-1,"Y"); |
259 |
gStyle->SetOptStat(0); |
260 |
|
261 |
// draw display area |
262 |
|
263 |
stringstream fromfile,title,hid,message; |
264 |
TString figsa="",figsav="",figsava=""; |
265 |
|
266 |
Int_t canvasx=1200; |
267 |
Int_t canvasy=900; |
268 |
figsav=out+ffile+"_FTrkScanQLook_EXPERT_ev"; |
269 |
figsav+=event+1; |
270 |
c1[e] = new TCanvas(figsav.Data(),"FTrkQLookSCAN",canvasx,canvasy); |
271 |
c1[e]->SetFillColor(10); |
272 |
c1[e]->Range(0,0,1,1); |
273 |
fromfile.str(""); |
274 |
fromfile<<"FTrkScanQLook_EXPERT File: "<<ffile<<" ----> Entry "<<event; |
275 |
t->DrawLatex(0.02,0.98,fromfile.str().c_str()); |
276 |
t->DrawLatex(0.60,0.98,cal.Data()); |
277 |
|
278 |
if(va1!=0){ |
279 |
figsava=out+ffile+"_FTrkScanQLook_EXPERT_ev"; |
280 |
figsava+=event+1; |
281 |
figsava+="_DSP"; |
282 |
figsava+=(int)(va1/100); |
283 |
figsava+="_VA1-"; |
284 |
figsava+=va1%100; |
285 |
cva[e] = new TCanvas(figsava.Data(),"TrkQLookSCAN VA1",canvasx,canvasy); |
286 |
cva[e]->SetFillColor(10); |
287 |
cva[e]->Range(0,0,1,1); |
288 |
fromfile.str(""); |
289 |
fromfile<<"FTrkScanQLook_EXPERT File: "<<ffile<<" ----> Entry "<<event<<" --> DSP "<<(int)(va1/100)<<" --> va1 "<<va1%100; |
290 |
t->DrawLatex(0.02,0.98,fromfile.str().c_str()); |
291 |
// t->DrawLatex(0.65,0.98,cal.Data()); |
292 |
} |
293 |
|
294 |
// draw pads |
295 |
TPad *trkpad[12],*pad; //pad for histos |
296 |
TPaveText *trkpadtext[12]; //pad for header |
297 |
Double_t posy = 0.95; // up y-coord - top pads |
298 |
Double_t hpad = 0.15; // pad height |
299 |
Double_t posx1=0; // left x-coord - pad column |
300 |
Double_t posx2=0; // right x-coord - pad olumn |
301 |
Double_t posx0=0; // x-coord - column division |
302 |
Double_t wrel = 0.6; // relative x size of first sub-column |
303 |
Double_t marg = 0.004; // margin among pads |
304 |
|
305 |
for(Int_t n = 0; n<12; n++){ |
306 |
if ( (n+1)%2 ) { |
307 |
if(n>1)posy = posy-(marg*2+hpad); |
308 |
posx1 = marg; |
309 |
posx2 = 0.5 - marg; |
310 |
posx0 = 0.5*wrel; |
311 |
} else { |
312 |
posx1 = posx1 + 0.5; |
313 |
posx2 = posx2 + 0.5; |
314 |
posx0 = posx0 + 0.5; |
315 |
}; |
316 |
|
317 |
/* -----------> pad for histograms */ |
318 |
trkpad[n] = new TPad("pad"," ",posx1,posy-hpad,posx0-marg,posy,18,0,0); |
319 |
trkpad[n]->SetFillColor(19); |
320 |
trkpad[n]->SetFrameFillColor(10); |
321 |
/* -----------> pad for header dump */ |
322 |
trkpadtext[n] = new TPaveText((posx0+marg),(posy-hpad),posx2,posy); |
323 |
/* -----------> HISTOGRAMS */ |
324 |
|
325 |
title<<"DSP "<<n+1; |
326 |
hid<<"h"<<n+e*100; |
327 |
histomax[n][e] = new TH1F(hid.str().c_str(),title.str().c_str(),3073,-0.5,3072.5); |
328 |
hid<<"hh"<<n+e*100; |
329 |
histocomp[n][e] = new TH1F(hid.str().c_str(),title.str().c_str(),3073,-0.5,3072.5); |
330 |
hid<<"hhh"<<n+e*100; |
331 |
histofull[n][e] = new TH1F(hid.str().c_str(),title.str().c_str(),3073,-0.5,3072.5); |
332 |
title.str(""); |
333 |
hid.str(""); |
334 |
} |
335 |
if(va1!=0){ |
336 |
title.str(""); |
337 |
title<<"DSP "<<(int)(va1/100)<<" -- va1 "<<va1%100; |
338 |
hid<<"va"<<e; |
339 |
histova[e] = new TH1F(hid.str().c_str(),title.str().c_str(),129,-0.5,128.5); |
340 |
} //end loop on views |
341 |
|
342 |
/* -----------> pad for histograms */ |
343 |
pad = new TPad("padva"," ",0,0,1,0.97,18,0,0); |
344 |
pad->SetFillColor(19); |
345 |
pad->SetFrameFillColor(10); |
346 |
|
347 |
// = = = = = = = = = = = = = = = = = = = = = = = = = |
348 |
// create header dump retrieving event info |
349 |
// = = = = = = = = = = = = = = = = = = = = = = = = = |
350 |
Int_t ndsp=0; |
351 |
|
352 |
Double_t whistomax[3072]; |
353 |
Double_t whisto[3072]; |
354 |
Double_t whistocomp[3072]; |
355 |
Double_t whistofull[3072]; |
356 |
|
357 |
//============================================= |
358 |
// transmitted words |
359 |
Int_t word = 0; |
360 |
Int_t iword = 0; |
361 |
Int_t TOTDATAlength_check = 0; |
362 |
Int_t ii=0,ifull[12],icomp[12],imax[12],nn=0; |
363 |
trkword thisword; |
364 |
|
365 |
Int_t address,ladder; |
366 |
|
367 |
for(Int_t n = 0; n<12; n++){ |
368 |
|
369 |
ndsp = trk->DSPnumber[n]; |
370 |
nn = ndsp-1; |
371 |
ifull[nn]=0; |
372 |
icomp[nn]=0; |
373 |
imax[nn]=0; |
374 |
if(ndsp>0){ |
375 |
if(ndsp<13){ |
376 |
|
377 |
/*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.* |
378 |
* |
379 |
* Write event LEVEL0 report |
380 |
* |
381 |
*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*/ |
382 |
|
383 |
trkpadtext[nn]->SetTextFont(40); |
384 |
trkpadtext[nn]->SetFillColor(33); |
385 |
trkpadtext[nn]->SetTextSize(0.012); |
386 |
trkpadtext[nn]->SetTextAlign(13); |
387 |
|
388 |
trkpadtext[nn]->AddText(" "); |
389 |
message<<"DAQ mode --------> "<<trk->DAQmode[n]; |
390 |
trkpadtext[nn]->AddText(message.str().c_str()); |
391 |
message.str(""); |
392 |
message<<"Event number --------> "<<trk->eventn[n]; |
393 |
trkpadtext[nn]->AddText(message.str().c_str()); |
394 |
message.str(""); |
395 |
message<<"13-bit words --------> "<<trk->DATAlength[n]; |
396 |
trkpadtext[nn]->AddText(message.str().c_str()); |
397 |
message.str(""); |
398 |
if (!(nn%2)&&trk->signcluster[n][0]!=0) message<<"L1 add: "<<trk->addrcluster[n][0]<<" - sign: "<<1024-(trk->signcluster[n][0]); |
399 |
else message<<"L1 add: "<<trk->addrcluster[n][0]<<" - sign: "<<(trk->signcluster[n][0]); |
400 |
trkpadtext[nn]->AddText(message.str().c_str()); |
401 |
message.str(""); |
402 |
if (!(nn%2)&&trk->signcluster[n][1]!=0) message<<"L2 add: "<<trk->addrcluster[n][1]<<" - sign: "<<1024-(trk->signcluster[n][1]); |
403 |
else message<<"L2 add: "<<trk->addrcluster[n][1]<<" - sign: "<<trk->signcluster[n][1]; |
404 |
trkpadtext[nn]->AddText(message.str().c_str()); |
405 |
message.str(""); |
406 |
if (!(nn%2)&&trk->signcluster[n][2]!=0) message<<"L3 add: "<<trk->addrcluster[n][2]<<" - sign: "<<1024-(trk->signcluster[n][2]); |
407 |
else message<<"L3 add: "<<trk->addrcluster[n][2]<<" - sign: "<<trk->signcluster[n][2]; |
408 |
trkpadtext[nn]->AddText(message.str().c_str()); |
409 |
message.str(""); |
410 |
message<<"NCLUST "<<trk->nclust[n]<<" CUTC "<<trk->cutc[n]<<" CUTCL "<<trk->cutcl[n]; |
411 |
trkpadtext[nn]->AddText(message.str().c_str()); |
412 |
message.str(""); |
413 |
message<<"Comp. time "<<trk->compressiontime[n]<<" x 0.051ms = "<<0.051*trk->compressiontime[n]<<" ms"; |
414 |
trkpadtext[nn]->AddText(message.str().c_str()); |
415 |
message.str(""); |
416 |
trkpadtext[nn]->AddText(" "); |
417 |
message<<"CRC --> "<<trk->crc[n]; |
418 |
trkpadtext[nn]->AddText(message.str().c_str()); |
419 |
message.str(""); |
420 |
message<<"FL1-6 --> "<<trk->fl1[n]<<" "<<trk->fl2[n]<<" "<<trk->fl3[n]<<" "<<trk->fl4[n]<<" "<<trk->fl5[n]<<" "<<trk->fl6[n]<<" FC "<<trk->fc[n]; |
421 |
trkpadtext[nn]->AddText(message.str().c_str()); |
422 |
message.str(""); |
423 |
trkpadtext[nn]->AddText(" "); |
424 |
trkpadtext[nn]->AddLine(0,0,0,0); |
425 |
message<<"PNum "<<trk->pnum[n]<<" - BId "<<trk->bid[n]<<" ALARM "<<trk->alarm[n]; |
426 |
trkpadtext[nn]->AddText(message.str().c_str()); |
427 |
message.str(""); |
428 |
message<<"Cmd "<<trk->cmdnum[n]<<" --- Answer length "<<trk->aswr[n]<<" byte"; |
429 |
trkpadtext[nn]->AddText(message.str().c_str()); |
430 |
message.str(""); |
431 |
trkpadtext[nn]->AddText(" "); |
432 |
|
433 |
TOTDATAlength_check = TOTDATAlength_check + trk->DATAlength[n]; |
434 |
|
435 |
/*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.* |
436 |
* |
437 |
* Plot event LEVEL0 histo |
438 |
* |
439 |
*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*.*/ |
440 |
|
441 |
//============================================= |
442 |
|
443 |
for(Int_t i = 0; i< 3072; i++){ |
444 |
whistomax[i] = -200; |
445 |
whistocomp[i] = -200; |
446 |
whistofull[i] = -200; |
447 |
whisto[i] = -200; |
448 |
} |
449 |
|
450 |
// =============== |
451 |
// trasmitted data |
452 |
// =============== |
453 |
|
454 |
address = 0; |
455 |
ladder = 1; |
456 |
for(Int_t i = 0; i < trk->DATAlength[n] ; i++){ |
457 |
word = trk->TrackerData.At(iword); |
458 |
thisword = datadecode(word); |
459 |
iword++; |
460 |
switch (thisword.type){ |
461 |
|
462 |
case 0: //ADC value |
463 |
whisto[address] = thisword.decode; |
464 |
address++; |
465 |
// cout << " adr " << address << "\n"; |
466 |
break; |
467 |
|
468 |
case 1: //address |
469 |
address = 1024*(ladder-1) + thisword.decode; |
470 |
// cout << " adr " << address << "\n"; |
471 |
break; |
472 |
|
473 |
case 2: //end-of-ladder |
474 |
ladder = thisword.decode; |
475 |
// cout << "Ladder " << ladder << "\n"; |
476 |
if(ladder==3){ |
477 |
// end of compressed data - FILL HISTO |
478 |
//cout << ">>> COMPRESSED data" << "\n"; |
479 |
for(ii = 0; ii < 3072; ii++){ |
480 |
whistocomp[ii]=whisto[ii]; |
481 |
whisto[ii] = -200; |
482 |
} |
483 |
address = 0; |
484 |
ladder = 1; |
485 |
}else if(ladder==6){ |
486 |
// end of full data - FILL HISTO |
487 |
//cout << ">>> FULL data" << "\n"; |
488 |
for(ii = 0; ii < 3072; ii++){ |
489 |
whistofull[ii]=whisto[ii]; |
490 |
whisto[ii] = -200; |
491 |
} |
492 |
address = 0; |
493 |
ladder = 1; |
494 |
}else{ |
495 |
if(ladder>3) ladder=ladder-3; |
496 |
address= ladder*1024; |
497 |
ladder = ladder + 1; |
498 |
} |
499 |
} |
500 |
} |
501 |
|
502 |
|
503 |
// =============== |
504 |
// maximum signals |
505 |
// =============== |
506 |
if(trk->signcluster[nn][0]!=0) whistomax[ 0+(int)trk->addrcluster[nn][0]] = whistocomp[ 0+(int)trk->addrcluster[nn][0]]; |
507 |
if(trk->signcluster[nn][1]!=0) whistomax[1024+(int)trk->addrcluster[nn][1]] = whistocomp[1024+(int)trk->addrcluster[nn][1]]; |
508 |
if(trk->signcluster[nn][2]!=0) whistomax[2048+(int)trk->addrcluster[nn][2]] = whistocomp[2048+(int)trk->addrcluster[nn][2]]; |
509 |
|
510 |
for(Int_t i = 0; i < 3072; i++){ |
511 |
if(whistomax[i]>-200) imax[nn]=1; |
512 |
if(whistocomp[i]>-200) icomp[nn]=1; |
513 |
if(whistofull[i]>-200) ifull[nn]=1; |
514 |
histomax[nn][e]->Fill((Float_t)i,whistomax[i]); |
515 |
histocomp[nn][e]->Fill((Float_t)i,whistocomp[i]); |
516 |
histofull[nn][e]->Fill((Float_t)i,whistofull[i]); |
517 |
if(va1!=0 && ndsp==(int)(va1/100) && i>=128*((va1%100)-1)){ |
518 |
if(i<128*(va1%100)) |
519 |
histova[e]->Fill((Float_t)(i-128*((va1%100)-1)+1),whistofull[i]); |
520 |
} |
521 |
} |
522 |
|
523 |
TBox b; |
524 |
TLine li,liva1; |
525 |
li.SetLineColor(38); |
526 |
li.SetLineStyle(4); |
527 |
li.SetLineWidth(2); |
528 |
liva1.SetLineColor(42); |
529 |
liva1.SetLineStyle(3); |
530 |
liva1.SetLineWidth(1); |
531 |
|
532 |
Float_t va1x=0; |
533 |
|
534 |
if(va1!=0){ |
535 |
cva[e]->cd(); |
536 |
pad->Draw(); |
537 |
pad->cd(); |
538 |
pad->SetFillColor(10); |
539 |
histova[e]->SetTitleSize(0.01); |
540 |
histova[e]->GetYaxis()->SetRangeUser(1500,4500); |
541 |
histova[e]->SetLineColor(40); |
542 |
histova[e]->SetFillColor(40); |
543 |
histova[e]->SetLineWidth(1); |
544 |
histova[e]->SetLineStyle(2); |
545 |
// histova[e]->GetYaxis()->SetLabelSize(0.03); |
546 |
histova[e]->Draw(""); |
547 |
cva[e]->Update(); |
548 |
} |
549 |
|
550 |
c1[e]->cd(); |
551 |
trkpadtext[nn]->Draw(); |
552 |
trkpad[nn]->Draw(); |
553 |
trkpad[nn]->cd(); |
554 |
trkpad[nn]->SetFillColor(10); |
555 |
|
556 |
histocomp[nn][e]->SetTitleSize(0.1); |
557 |
histocomp[nn][e]->GetYaxis()->SetRangeUser(-500,4500); |
558 |
histocomp[nn][e]->SetLineStyle(1); |
559 |
histocomp[nn][e]->SetLineColor(38); |
560 |
histocomp[nn][e]->SetFillColor(38); |
561 |
histocomp[nn][e]->SetLineWidth(1); |
562 |
histocomp[nn][e]->Draw(""); |
563 |
|
564 |
histofull[nn][e]->SetLineColor(40); |
565 |
histofull[nn][e]->SetFillColor(40); |
566 |
histofull[nn][e]->SetLineWidth(1); |
567 |
histofull[nn][e]->SetLineStyle(2); |
568 |
|
569 |
histomax[nn][e]->SetLineColor(2); |
570 |
histomax[nn][e]->SetLineWidth(1); |
571 |
histomax[nn][e]->SetLineStyle(3); |
572 |
|
573 |
if(ifull[nn]==1) histofull[nn][e]->Draw("9bsame]["); |
574 |
if(icomp[nn]==1) histocomp[nn][e]->Draw("9bsame]["); |
575 |
if(imax[nn]==1) histomax[nn][e]->Draw("same]["); |
576 |
histocomp[nn][e]->Draw("axis same"); |
577 |
if(nn==0){ |
578 |
b.SetFillColor(107); |
579 |
b.SetFillStyle(3945); |
580 |
b.DrawBox(768.,-500.,2047.,4500.); |
581 |
} |
582 |
else if(nn==1){ |
583 |
b.SetFillColor(6); |
584 |
b.SetFillStyle(3945); |
585 |
b.DrawBox(2944.,-500.,3060.,4500.); |
586 |
|
587 |
b.SetFillColor(107); |
588 |
b.SetFillStyle(3954); |
589 |
//b.DrawBox(384.,-500.,512.,4500.); |
590 |
b.DrawBox(2816.,-500.,2944.,4500.); |
591 |
b.DrawBox(2048.,-500.,2176.,4500.); |
592 |
} |
593 |
else if(nn==4){ |
594 |
b.SetFillColor(107); |
595 |
b.SetFillStyle(3954); |
596 |
b.DrawBox(384.,-500.,512.,4500.); |
597 |
} |
598 |
else if(nn==6){ |
599 |
b.SetFillColor(6); |
600 |
b.SetFillStyle(3945); |
601 |
b.DrawBox(2560.,-500.,2816.,4500.); |
602 |
b.DrawBox(1024.,-500.,1535.,4500.); |
603 |
|
604 |
b.SetFillColor(107); |
605 |
b.SetFillStyle(3954); |
606 |
b.DrawBox(512.,-500.,768.,4500.); |
607 |
b.DrawBox(1536.,-500.,1792.,4500.); |
608 |
} |
609 |
else if(nn==7){ |
610 |
b.SetFillColor(107); |
611 |
b.SetFillStyle(3954); |
612 |
b.DrawBox(512.,-500.,768.,4500.); |
613 |
} |
614 |
else if(nn==8){ |
615 |
b.SetFillColor(107); |
616 |
b.SetFillStyle(3954); |
617 |
b.DrawBox(512.,-500.,768.,4500.); |
618 |
} |
619 |
else if(nn==9){ |
620 |
b.SetFillColor(107); |
621 |
b.SetFillStyle(3954); |
622 |
b.DrawBox(256.,-500.,384.,4500.); |
623 |
//b.DrawBox(1280.,-500.,1408.,4500.); |
624 |
//b.DrawBox(1792.,-500.,1920.,4500.); |
625 |
} |
626 |
else if(nn==10){ |
627 |
b.SetFillColor(107); |
628 |
b.SetFillStyle(3954); |
629 |
b.DrawBox(2048.,-500.,3070.,4500.); |
630 |
} |
631 |
else if(nn==11){ |
632 |
b.SetFillColor(6); |
633 |
b.SetFillStyle(3945); |
634 |
b.DrawBox(768.,-500.,1024.,4500.); |
635 |
|
636 |
b.SetFillColor(107); |
637 |
b.SetFillStyle(3954); |
638 |
b.DrawBox(0.,-500.,512.,4500.); |
639 |
b.DrawBox(1920.,-500.,2560.,4500.); |
640 |
} |
641 |
for(int va=1; va<24; va++){ |
642 |
va1x=128*va; |
643 |
liva1.DrawLine(va1x,-500.,va1x,4500.); |
644 |
} |
645 |
li.DrawLine(1024.5,-500.,1024.5,4500.); |
646 |
li.DrawLine(2048.5,-500.,2048.5,4500.); |
647 |
c1[e]->Update(); |
648 |
} |
649 |
} |
650 |
|
651 |
}//end loop on views |
652 |
|
653 |
stringstream nom1,nom2,nom3; |
654 |
nom1<<out<<ffile<<"_FTrkScanQLook_EXPERT.ps("; |
655 |
nom2<<out<<ffile<<"_FTrkScanQLook_EXPERT.ps"; |
656 |
nom3<<out<<ffile<<"_FTrkScanQLook_EXPERT.ps)"; |
657 |
|
658 |
if(!strcmp(outfile.Data(),"ps")||!strcmp(outfile.Data(),"pdf")){ |
659 |
if(e==0){ |
660 |
c1[e]->Print(nom1.str().c_str(),"Portrait"); |
661 |
if(va1!=0) cva[e]->Print(nom2.str().c_str(),"Portrait"); |
662 |
} |
663 |
if(e>0 && tot>2){ |
664 |
c1[e]->Print(nom2.str().c_str(),"Portrait"); |
665 |
if(va1!=0) cva[e]->Print(nom2.str().c_str(),"Portrait"); |
666 |
} |
667 |
if(e==tot-1){ |
668 |
c1[e]->Print(nom2.str().c_str(),"Portrait"); |
669 |
if(va1!=0) cva[e]->Print(nom3.str().c_str(),"Portrait"); |
670 |
} |
671 |
} |
672 |
else{ |
673 |
figsav+="."+outfile; |
674 |
c1[e]->Print(figsav.Data()); |
675 |
if(va1!=0){ |
676 |
figsava+="."+outfile; |
677 |
cva[e]->Print(figsava.Data()); |
678 |
} |
679 |
} |
680 |
} |
681 |
|
682 |
// |
683 |
// Convert ps to pdf if required |
684 |
if(!strcmp(outfile.Data(),"pdf")){ |
685 |
stringstream com; |
686 |
com<<"ps2pdf13 "<<out<<ffile<<"_FTrkScanQLook_EXPERT.ps "<<out<<ffile<<"_FTrkScanQLook_EXPERT.pdf"; |
687 |
system(com.str().c_str()); |
688 |
printf("\n---> ps file converted in pdf format!\n"); |
689 |
com.str(""); |
690 |
com<<"rm -f "<<out<<ffile<<"_FTrkScanQLook_EXPERT.ps "; |
691 |
system(com.str().c_str()); |
692 |
printf("---> ps file removed!\n\n"); |
693 |
com.str(""); |
694 |
} |
695 |
|
696 |
return; |
697 |
} |