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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.19 - (show annotations) (download)
Fri Jan 17 15:10:40 2014 UTC (10 years, 10 months ago) by mocchiut
Branch: MAIN
CVS Tags: HEAD
Changes since 1.18: +4 -0 lines
Error occurred while calculating annotation data.
Compilation warnings using GCC4.7 fixed

1 /**
2 * FTrkQLook_BASIC.cxx
3 *
4 * autor: D.Fedele
5 * version v1r12
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 <TPaveText.h>
16 #include <TLatex.h>
17 #include <TCanvas.h>
18 #include <TGraph.h>
19 #include <TStyle.h>
20 #include <TFile.h>
21 #include <TTree.h>
22 #include <TArrow.h>
23 #include <TAxis.h>
24 #include <TROOT.h>
25 #include <cstdlib>
26 #include <sys/stat.h>
27 //
28 #include <physics/tracker/TrackerEvent.h>
29 #include <PscuHeader.h>
30 #include <EventHeader.h>
31 #include <RunHeaderEvent.h>
32 #include <EventCounter.h>
33 #include <PacketType.h>
34 //
35 #define MAXSTORAGE 50000
36
37 void FTrkQLook_BASIC(TString file,Int_t fromevent,Int_t toevent, TString outdir,TString outfile)
38 {
39 //
40 // obtain information about the data file and select the output dir
41 Int_t dwpos = file.Last('/');
42 Int_t dwpos1 = file.Last('.');
43 TString base,ffile ;
44 ffile=file(dwpos+1,dwpos1-(dwpos+1));
45 if(dwpos>0) base=file(0,dwpos);
46
47 TString out;
48 if(outdir.Length()==0){
49 out = base;
50 }else{
51 out = outdir;
52 }
53 if(out.Last('/')+1<out.Length()) out+="/";
54
55 //
56 // inizialise the variables and open the file
57 pamela::tracker::TrackerEvent *te=0;
58 pamela::EventHeader *eh=0,*eH=0,*ceh=0;
59 pamela::RunHeaderEvent *reh=0;
60 pamela::PscuHeader *ph=0,*pH=0;
61 pamela::EventCounter *cod=0;
62
63 pamela::PacketType *pctp=0;
64
65 TFile *datafile = new TFile(file);
66 TTree *otr = (TTree*)datafile->Get("RunHeader");
67 otr->SetBranchAddress("Header",&eH);
68 otr->SetBranchAddress("RunHeader",&reh);
69 TTree *tr = (TTree*)datafile->Get("Physics");
70 tr->SetBranchAddress("Tracker",&te);
71 tr->SetBranchAddress("Header",&eh);
72 TTree *ctr = (TTree*)datafile->Get("CalibTrk1");
73 ctr->SetBranchAddress("Header",&ceh);
74
75 Long64_t neventC = ctr->GetEntries();
76 Long64_t nevent = tr->GetEntries();
77 Long64_t neventH = otr->GetEntries();
78 Int_t minevent=0;
79 Int_t maxevent=0;
80
81 printf("Number of total events: %lld\nNumber of total header events: %lld\n",nevent,neventH);
82 printf("Number of calibration events: %lld\n",neventC);
83
84 if (nevent<=0){
85 datafile->Close();
86 return;
87 }
88 if ( fromevent > toevent && toevent>0 ){
89 printf("It must be fromevent < toevent \n");
90 return;
91 }
92 if ( fromevent > nevent || fromevent < 0 ) {
93 printf("You can choose fromevent between 0 (all) and %lld \n",nevent);
94 return;
95 }
96 if ( toevent > nevent || toevent < 0 ) {
97 printf("You can choose toevent between 0 (all) and %lld \n",nevent);
98 return;
99 }
100 if ( fromevent == 0 ) {
101 minevent = 0;
102 maxevent = nevent;
103 } else {
104 minevent = fromevent;
105 if ( toevent > 0 ){
106 maxevent = toevent+1;
107 } else if (toevent > nevent) {
108 maxevent = nevent;
109 } else {
110 maxevent = nevent;
111 }
112 nevent=maxevent-minevent ;
113 }
114
115 //
116 // other variables definitions
117 stringstream oss,fromfile,isfile;
118 //
119 // information about the RunHeader
120 ULong_t HOBT[neventH];
121 Int_t trk_cal_us[neventH];
122 for (Int_t vi=0; vi<neventH;vi++){
123 HOBT[vi]=0;
124 trk_cal_us[vi]=0;
125 }
126
127 //
128 // information about RunHeader
129 Int_t countnboot=1;
130 for (Int_t ev=0; ev<neventH; ev++){
131 otr->GetEntry(ev);
132 pH = eH->GetPscuHeader();
133 HOBT[ev]= pH->GetOrbitalTime();
134 trk_cal_us[ev]=reh->TRK_CALIB_USED;
135 if((HOBT[ev]<HOBT[ev-1]) && ev>0)
136 countnboot+=1;
137 }
138 countnboot+=(Int_t)(2*nevent/MAXSTORAGE);
139 // printf("\ncountnboot=%d\n",countnboot);
140
141 //
142 // information about calibration OBT
143 ULong_t COBT[neventC];
144 for (Int_t vi=0; vi<neventC;vi++){
145 COBT[vi]=0;
146 }
147 for (Int_t ev=0; ev<neventC; ev++){
148 ctr->GetEntry(ev);
149 pH = ceh->GetPscuHeader();
150 COBT[ev]= pH->GetOrbitalTime();
151 }
152
153 //
154 // Style options
155 gStyle->SetLabelSize(0.06,"x");
156 gStyle->SetLabelSize(0.06,"y");
157 gStyle->SetStatFontSize(0.075);
158 gStyle->SetOptStat(10);
159 gStyle->SetFillColor(10);
160 gStyle->SetTitleFontSize(0.1);
161 gStyle->SetTitleFillColor(10);
162 gStyle->SetTitleOffset(0.8,"y");
163 gStyle->SetTitleOffset(0.9,"x");
164 gStyle->SetTitleSize(0.06,"y");
165 gStyle->SetTitleSize(0.055,"x");
166
167 //***************************************************************************************
168 // LOOP on each event
169 //***************************************************************************************
170
171 if (fromevent!=0)
172 printf("\n Scan of events from %i to %i ... \n",minevent,maxevent-1);
173 else
174 printf("\n Scan of events from %i to %i ... \n",minevent+1,maxevent);
175
176
177 Int_t minev=minevent,maxev=maxevent,hin=0,hfin=0,cin=0,cfin=0;
178 TPad *pad[12][countnboot];
179 TGraph *dataletime[12][countnboot],*dataletime1[12][countnboot];
180 TCanvas *DataTimeCanv[countnboot];
181 for(Int_t ii=0; ii<countnboot;ii++){
182 fromfile<<"FTrkQLook_BASIC File: "<<ffile;
183 isfile<<"DATALENGTH vs. OBT pag"<<ii+1;
184 DataTimeCanv[ii]=new TCanvas(isfile.str().c_str(),isfile.str().c_str(),900,1200);
185 DataTimeCanv[ii]->SetFillColor(10);
186 DataTimeCanv[ii]->Range(0,0,100,100);
187
188 TLatex *t=new TLatex();
189 t->SetTextFont(32);
190 t->SetTextColor(1);
191 t->SetTextAlign(12);
192 t->SetTextSize(0.02);
193 t->DrawLatex(2.,98.7,fromfile.str().c_str());
194 TLatex *t1=new TLatex();
195 t1->SetTextFont(32);
196 t1->SetTextColor(1);
197 t1->SetTextAlign(12);
198 t1->SetTextSize(0.02);
199 t1->DrawLatex(70.,98.7,isfile.str().c_str());
200 fromfile.str("");
201 isfile.str("");
202
203 fromfile<<"D = Default Calibration";
204 isfile<<"O = OnLine Calibration";
205 t->SetTextColor(6);
206 t->SetTextSize(0.018);
207 t->DrawLatex(70.,97.,fromfile.str().c_str());
208 t->SetTextColor(3);
209 t->DrawLatex(70.,96.,isfile.str().c_str());
210 fromfile.str("");
211 isfile.str("");
212
213 fromfile<<"The green arrow (if present) points out the time of the online calibration";
214 t->DrawLatex(7.,96.,fromfile.str().c_str());
215 fromfile.str("");
216
217 //*************************************************************************************
218 //book pads and histos
219 //***************************************************************************************
220
221
222 Float_t posy = 0.95; // up y-coord - top pads
223 Float_t hpad = 0; // pad height
224 Float_t posx1=0; // left x-coord - pad column
225 Float_t posx0=0; // x-coord - column division
226 Float_t wrel = 0; // relative x size of first sub-column
227 Float_t marg = 0.004; // margin among pads
228
229 hpad = (posy-marg*11)/6;
230 wrel = (1-marg*4)/2;
231 stringstream title;
232 stringstream hid;
233
234 for(Int_t n = 0; n<12; n++) {
235 if ( (n+1)%2==1 ) {
236 if(n>1) posy = posy-(marg*2+hpad);
237 posx1 = marg;
238 posx0 = posx1 + wrel;
239 }
240 else {
241 posx1 = posx0 + 2*marg;
242 posx0 = posx1 + wrel;
243 }
244
245 /* -----------> pad for histograms */
246 oss<<"pad"<<n*100+ii;
247 pad[n][ii]=new TPad(oss.str().c_str()," ",posx1,posy-hpad,posx0,posy,18,0,0);
248 oss.str("");
249 };
250
251 TLine li;
252 li.SetLineColor(1);
253 li.SetLineStyle(1);
254 li.SetLineWidth(1);
255
256 TArrow ar;
257 ar.SetLineColor(3);
258 stringstream calus;
259
260 TLatex *t2=new TLatex();
261 t2->SetTextFont(32);
262 t2->SetTextColor(1);
263 t2->SetTextAlign(13);
264 t2->SetTextSize(0.08);
265
266 Int_t i=0;
267 ULong_t x[MAXSTORAGE], xb[MAXSTORAGE];
268 Float_t yyd[MAXSTORAGE][12],yyb[MAXSTORAGE][12];
269 Int_t countbad[12];
270 Float_t perc=0;
271 Double_t xMIN=0.,xMAX=0.;
272 for (Int_t n=0; n<12 ; n++)
273 countbad[n]=0;
274
275 //
276 // obtain values of the datalenght
277 for (Int_t ev=minev; ev<maxevent; ev++){
278 tr->GetEntry(ev);
279 ph = eh->GetPscuHeader();
280 cod = eh->GetCounter();
281
282 if(ev==minev){
283 if(cod->Get(pctp->CalibTrk1)>0) cin=cod->Get(pctp->CalibTrk1)-1;
284 else cin=cod->Get(pctp->CalibTrk1);
285 if(cin==cfin-1) cin+=1;
286
287 if(cod->Get(pctp->RunHeader)>0) hin=cod->Get(pctp->RunHeader)-1;
288 else hin=cod->Get(pctp->RunHeader);
289 if(hin==hfin-1) hin+=1;
290 }
291
292 if(ev==maxevent-1) maxev=maxevent-1;
293
294 if(((ULong_t)ph->GetOrbitalTime()<x[ev-minev-1] && ev-minev!=0) || ev-minev==MAXSTORAGE){
295 maxev=ev;
296 break;
297 }
298 else{
299 cfin=cod->Get(pctp->CalibTrk1);
300 hfin=cod->Get(pctp->RunHeader);
301 x[(ev-minev)]= ph->GetOrbitalTime();
302 i=0;
303
304 for (Int_t n=0; n<12 ; n++){
305 i=te->DSPnumber[n]-1;
306 if(i>=0){
307 if(i<12){
308 yyb[countbad[i]][i]=0;
309 xb[countbad[i]]= 0;
310
311 yyd[(ev-minev)][i]=te->DATAlength[n];
312 if(i==6){
313 if(yyd[(ev-minev)][i]>1500){
314 if(yyd[(ev-minev)][i]<3075){
315 yyb[countbad[i]][i]= yyd[(ev-minev)][i];
316 xb[countbad[i]]= x[(ev-minev)];
317 countbad[i]+=1;
318 }
319 }
320 }
321 else{
322 if(yyd[(ev-minev)][i]>750){
323 if(yyd[(ev-minev)][i]<3075){
324 yyb[countbad[i]][i]= yyd[(ev-minev)][i];
325 xb[countbad[i]]= x[(ev-minev)];
326 countbad[i]+=1;
327 }
328 }
329 }
330 }
331 }
332 }
333 }
334 }
335
336 //
337 // define limit for the Xaxis of the graphs
338
339 xMAX=x[maxev-minev-1]+(x[maxev-minev-1]-x[0])/10;
340 if(xMAX>1000000) xMIN=x[0]-(x[maxev-minev-1]-x[0])/10;
341 if(xMAX<1000000 || xMIN<0) xMIN=0.;
342 if(xMIN==xMAX) xMIN=x[0] * x[0]/xMAX;
343
344
345 //
346 // Draw Histos
347 for (Int_t n=0; n<12 ; n++){
348
349 Double_t yd[maxev-minev],ybd[maxev-minev];
350 Double_t xd[maxev-minev],xbd[maxev-minev];
351
352 for(Int_t v=0;v<maxev-minev;v++){
353 yd[v]=(Double_t)yyd[v][n];
354 ybd[v]=(Double_t)yyb[v][n];
355 xd[v]=(Double_t)x[v];
356 xbd[v]=(Double_t)xb[v];
357 }
358
359 if(xMAX<1000000){
360 xd[maxev-minev-1]=0;
361 yd[maxev-minev-1]=0.;
362 }
363
364 if((maxev-minev)>1000){
365 perc=(countbad[n]*100)/(maxev-minev);
366 if(perc>10) pad[n][ii]->SetFillColor(2);
367 else pad[n][ii]->SetFillColor(10);
368
369 }
370 else{
371 if(countbad[n]>=100) pad[n][ii]->SetFillColor(2);
372 else pad[n][ii]->SetFillColor(10);
373 }
374
375 oss<<"DSP "<<n+1;
376 DataTimeCanv[ii]->cd();
377 pad[n][ii]->SetFrameFillColor(10);
378 pad[n][ii]->Draw();
379 pad[n][ii]->cd();
380 dataletime[n][ii]= new TGraph((maxev-minev),xd,yd);
381 if(maxev-minev==1) dataletime[n][ii]->GetXaxis()->SetNdivisions(-501);
382 dataletime[n][ii]->SetTitle(oss.str().c_str());
383 dataletime[n][ii]->GetXaxis()->SetTitle("OBT (ms)");
384 dataletime[n][ii]->GetXaxis()->CenterTitle();
385 dataletime[n][ii]->GetXaxis()->SetRangeUser(xMIN,xMAX);
386 dataletime[n][ii]->GetYaxis()->SetTitle("datalength (Word 13 bit)");
387 dataletime[n][ii]->GetYaxis()->CenterTitle();
388 if(i==6) dataletime[n][ii]->GetYaxis()->SetRangeUser(0,4500);
389 else dataletime[n][ii]->GetYaxis()->SetRangeUser(0,4000);
390 dataletime[n][ii]->SetMarkerStyle(21);
391 if((maxev-minev)<50) dataletime[n][ii]->SetMarkerSize(0.5);
392 else dataletime[n][ii]->SetMarkerSize(0.1);
393 dataletime[n][ii]->SetMarkerColor(4);
394 dataletime[n][ii]->Draw("ap");
395
396
397 if((maxev-minev)>1000 && perc>10){
398 dataletime1[n][ii]= new TGraph(countbad[n],xbd,ybd);
399 dataletime1[n][ii]->SetMarkerStyle(21);
400 dataletime1[n][ii]->SetMarkerSize(0.1);
401 dataletime1[n][ii]->SetMarkerColor(2);
402 dataletime1[n][ii]->Draw("psame");
403 }
404 else if((maxev-minev)<1000 && countbad[n]>=100){
405 dataletime1[n][ii]= new TGraph(countbad[n],xbd,ybd);
406 dataletime1[n][ii]->SetMarkerStyle(21);
407 if((maxev-minev)<50) dataletime1[n][ii]->SetMarkerSize(0.5);
408 else dataletime1[n][ii]->SetMarkerSize(0.1);
409 dataletime1[n][ii]->SetMarkerColor(2);
410 dataletime1[n][ii]->Draw("psame");
411 }
412 li.SetLineColor(1);
413 li.SetLineStyle(1);
414 li.SetLineWidth(1);
415 if(i!=6) li.DrawLine(xMIN,750,xMAX,750);
416 else li.DrawLine(xMIN,1500,xMAX,1500);
417 li.DrawLine(xMIN,3075,xMAX,3075);
418
419 li.SetLineColor(12);
420 li.SetLineStyle(4);
421 li.SetLineWidth(1);
422 for(Int_t j=hin;j<hfin;j++){
423 if(i==6) li.DrawLine(HOBT[j],0.,HOBT[j],4500.);
424 else li.DrawLine(HOBT[j],0.,HOBT[j],4000.);
425 if(trk_cal_us[j]==104){
426 calus<<"D";
427 t2->SetTextColor(6);
428 if(i==6) t2->DrawLatex(HOBT[j],4350.,calus.str().c_str());
429 else t2->DrawLatex(HOBT[j],3850.,calus.str().c_str());
430 calus.str("");
431 }
432 else{
433 calus<<"O";
434 t2->SetTextColor(3);
435 if(i==6) t2->DrawLatex(HOBT[j],4350.,calus.str().c_str());
436 else t2->DrawLatex(HOBT[j],3850.,calus.str().c_str());
437 calus.str("");
438 }
439 }
440 for(Int_t j=cin;j<cfin;j++){
441 if(i==6) ar.DrawArrow(COBT[j],1700.,COBT[j],2700.,0.01,"<");
442 else ar.DrawArrow(COBT[j],1000.,COBT[j],2000.,0.01,"<");
443 }
444
445 oss.str("");
446 DataTimeCanv[ii]->Update();
447 }
448
449 minev=maxev;
450 if(maxev==maxevent-1) {
451 countnboot=ii+1;
452 break;
453 }
454 }
455 printf("... end of packets. \n");
456
457 //*************************************************************************
458 // Save output Files
459 //*************************************************************************
460 stringstream nom1,nom2,nom3;
461
462 for(Int_t fl=0;fl<countnboot;fl++){
463 if(countnboot==1){
464 nom1<<ffile<<"_FTrkQLook_BASIC."<<outfile.Data();
465 DataTimeCanv[fl]->Print(out+nom1.str().c_str());
466 nom1.str("");
467 }
468
469 if(countnboot>=2){
470 if(!strcmp(outfile.Data(),"ps") || !strcmp(outfile.Data(),"pdf")){
471 nom1.str("");
472 nom2.str("");
473 nom3.str("");
474 nom1<<ffile<<"_FTrkQLook_BASIC.ps(";
475 nom2<<ffile<<"_FTrkQLook_BASIC.ps";
476 nom3<<ffile<<"_FTrkQLook_BASIC.ps)";
477 if(fl==0) DataTimeCanv[fl]->Print(out+nom1.str().c_str(),"portrait");
478 else if(fl==countnboot-1) DataTimeCanv[fl]->Print(out+nom3.str().c_str(),"portrait");
479 else DataTimeCanv[fl]->Print(out+nom2.str().c_str(),"portrait");
480
481 }
482 else{
483 nom1.str("");
484 nom1<<ffile<<"_FTrkQLook_BASIC-pag"<<fl+1<<"."<<outfile.Data();
485 DataTimeCanv[fl]->Print(out+nom1.str().c_str());
486 }
487 }
488 }
489
490 //
491 // Convert ps to pdf if required
492 if(!strcmp(outfile.Data(),"pdf") && countnboot>=2){
493 stringstream com;
494 com<<"ps2pdf13 "<<out<<ffile<<"_FTrkQLook_BASIC.ps "<<out<<ffile<<"_FTrkQLook_BASIC.pdf";
495 system(com.str().c_str());
496 printf("\n---> ps file converted in pdf format!\n");
497 com.str("");
498 com<<"rm -f "<<out<<ffile<<"_FTrkQLook_BASIC.ps";
499 system(com.str().c_str());
500 printf("---> ps file removed!\n\n");
501 com.str("");
502 }
503
504 datafile->Close();
505 gROOT->Reset();
506 return;
507 }

  ViewVC Help
Powered by ViewVC 1.1.23