/[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.2 - (show annotations) (download)
Thu Jun 1 15:23:12 2006 UTC (18 years, 6 months ago) by pam-fi
Branch: MAIN
Changes since 1.1: +345 -184 lines
*** empty log message ***

1 /**
2 * FTrkQLook_BASIC.cpp
3 *
4 * autor: D.Fedele
5 * version 2.0
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 <TTree.h>
21 #include <TArrow.h>
22 //
23 #include <physics/tracker/TrackerEvent.h>
24 #include <PscuHeader.h>
25 #include <EventHeader.h>
26 #include <RunHeaderEvent.h>
27 //
28
29 void stringcopy(TString& s1, const TString& s2, Int_t from=0, Int_t to=0){
30 if ( to == 0 ){
31 Int_t t2length = s2.Length();
32 s1 = "";
33 to = t2length;
34 };
35 for (Int_t i = from; i<to; i++){
36 s1.Append(s2[i],1);
37 };
38 };
39
40
41 void FTrkQLook_BASIC(TString file,Int_t fromevent,Int_t toevent, TString outdir,TString outfile)
42 {
43 //
44 // obtain information about the data file and select the output dir
45 const string filepath=file.Data();
46 Int_t dwpos = filepath.rfind("/");
47 Int_t dwpos1 = filepath.find(".root");
48 TString fpath=(filepath.c_str());
49 TString base,ffile ;
50 stringcopy(ffile,fpath,dwpos+1,dwpos1);
51 stringcopy(base,fpath,0,dwpos);
52 if(dwpos>0) base+="/";
53
54 TString out;
55 if(outdir.Length()==0){
56 out = base;
57 }else{
58 out = outdir;
59 }
60
61 //
62 // inizialise the variables and open the file
63 pamela::tracker::TrackerEvent *te=0;
64 pamela::EventHeader *eh=0,*eH=0,*ceh=0;
65 pamela::RunHeaderEvent *reh=0;
66 pamela::PscuHeader *ph=0,*pH=0;
67
68 TFile *datafile = new TFile(file);
69 TTree *otr = (TTree*)datafile->Get("RunHeader");
70 otr->SetBranchAddress("Header",&eH);
71 otr->SetBranchAddress("RunHeader",&reh);
72 TTree *tr = (TTree*)datafile->Get("Physics");
73 tr->SetBranchAddress("Tracker",&te);
74 tr->SetBranchAddress("Header",&eh);
75 TTree *ctr = (TTree*)datafile->Get("CalibTrk1");
76 ctr->SetBranchAddress("Header",&ceh);
77
78
79 Long64_t neventC = ctr->GetEntries();
80 Long64_t nevent = tr->GetEntries();
81 Long64_t neventH = otr->GetEntries();
82 Int_t minevent=0;
83 Int_t maxevent=0;
84
85 printf("Number of total events: %lld\nNumber of total header events: %lld\n",nevent,neventH);
86 printf("Number of calibration events: %lld\n",neventC);
87
88 if (nevent<=0){
89 datafile->Close();
90 return;
91 }
92 if ( fromevent > toevent && toevent>0 ){
93 printf("It must be fromevent < toevent \n");
94 return;
95 }
96 if ( fromevent > nevent || fromevent < 0 ) {
97 printf("You can choose fromevent between 0 (all) and %lld \n",nevent);
98 return;
99 }
100 if ( toevent > nevent || toevent < 0 ) {
101 printf("You can choose toevent between 0 (all) and %lld \n",nevent);
102 return;
103 }
104 if ( fromevent == 0 ) {
105 minevent = 0;
106 maxevent = nevent;
107 } else {
108 minevent = fromevent;
109 if ( toevent > 0 ){
110 maxevent = toevent+1;
111 } else if (toevent > nevent) {
112 maxevent = nevent;
113 } else {
114 maxevent = toevent+1;
115 }
116 nevent=maxevent-minevent ;
117 }
118
119 //
120 // other variables definitions
121 stringstream oss,fromfile,isfile;
122 const Int_t sizeH=neventH;
123 const Int_t sizeC=neventC;
124 Int_t count=0,trk_cal_us[sizeH],countrun=1;
125 Double_t perc=0,xMIN=0,xMAX=0;
126 ULong64_t HOBT[sizeH],COBT[sizeC];
127
128 for (Int_t vi=0; vi<sizeH;vi++){
129 HOBT[vi]=0;
130 trk_cal_us[vi]=0;
131 }
132 for (Int_t vi=0; vi<sizeC;vi++){
133 COBT[vi]=0;
134 }
135
136 //***************************************************************************************
137 // LOOP on each event
138 //***************************************************************************************
139
140 //
141 // information about trk_calib_used
142 for (Int_t ev=0; ev<neventH; ev++){
143 otr->GetEntry(ev);
144 pH = eH->GetPscuHeader();
145 HOBT[ev]= pH->GetOrbitalTime();
146 trk_cal_us[ev]=reh->TRK_CALIB_USED;
147 if((HOBT[ev]<HOBT[ev-1]) && ev>0)
148 countrun+=1;
149 // printf("\n%lld\t\tcountrun=%d\n",HOBT[ev],countrun);
150 }
151 countrun+=20;
152 printf("\ncountrun=%d\n",countrun);
153 //
154 // information about calibration OBT
155 for (Int_t ev=0; ev<neventC; ev++){
156 ctr->GetEntry(ev);
157 pH = ceh->GetPscuHeader();
158 COBT[ev]= pH->GetOrbitalTime();
159 }
160
161 //****************************************************************************************
162 //Output figure
163 //****************************************************************************************
164 gStyle->SetLabelSize(0.06,"x");
165 gStyle->SetLabelSize(0.06,"y");
166 gStyle->SetStatFontSize(0.075);
167 gStyle->SetOptStat(10);
168 gStyle->SetFillColor(10);
169 gStyle->SetTitleFontSize(0.1);
170 gStyle->SetTitleFillColor(10);
171 gStyle->SetTitleOffset(0.8,"y");
172 gStyle->SetTitleOffset(0.9,"x");
173 gStyle->SetTitleSize(0.06,"y");
174 gStyle->SetTitleSize(0.055,"x");
175
176
177 Int_t minev=minevent,maxev=maxevent,countHOBT=0,countCOBT=0;
178 TPad *pad[12][countrun] ; //pad for histos
179 TGraph *dataletime[12][countrun],*dataletime1[12][countrun];
180 TCanvas *DataTimeCanv[countrun];
181 for(Int_t ii=0; ii<countrun;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(10.,96.,fromfile.str().c_str());
215 fromfile.str("");
216
217 //*************************************************************************************
218 //book pads and histos
219 //***************************************************************************************
220
221
222 Double_t posy = 0.95; // up y-coord - top pads
223 Double_t hpad = 0; // pad height
224 Double_t posx1=0; // left x-coord - pad column
225 Double_t posx0=0; // x-coord - column division
226 Double_t wrel = 0; // relative x size of first sub-column
227 Double_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 //**********************************************************************************
259 // Fill Graphs and Histos
260 //**********************************************************************************
261 stringstream calus;
262
263 TLatex *t2=new TLatex();
264 t2->SetTextFont(32);
265 t2->SetTextColor(1);
266 t2->SetTextAlign(13);
267 t2->SetTextSize(0.08);
268
269 Int_t i=0;
270 Double_t x[30000], xb[30000];
271 Double_t yyd[30000][12],yyb[30000][12];
272 for (Int_t ev=minev; ev<maxevent; ev++){
273 tr->GetEntry(ev);
274 ph = eh->GetPscuHeader();
275
276 if(ev==maxevent-1) maxev=maxevent-1;
277 // printf("prova00\n");
278 if((ph->GetOrbitalTime()<x[ev-minev-1]&&ev-minev!=0) || ev-minev>30000){
279 maxev=ev;
280 // printf("prova12\n");
281 break;
282 }
283 else{
284 x[(ev-minev)]= ph->GetOrbitalTime();
285 i=0;
286
287 // printf("prova0 %d\n ",ev);
288 for (Int_t n=0; n<12 ; n++){
289 perc=0;
290 count=0;
291 yyb[count][i]=0;
292 xb[count]= 0;
293
294 i=te->DSPnumber[n]-1;
295
296 yyd[(ev-minev)][i]=te->DATAlength[n];
297 if(i==6){
298 if(yyd[(ev-minev)][i]>1500){
299 if(yyd[(ev-minev)][i]<3075){
300 yyb[count][i]= yyd[(ev-minev)][i];
301 xb[count]= x[(ev-minev)];
302 count++;
303 }
304 }
305 }
306 else{
307 if(yyd[(ev-minev)][i]>750){
308 if(yyd[(ev-minev)][i]<3075){
309 yyb[count][i]= yyd[(ev-minev)][i];
310 xb[count]= x[(ev-minev)];
311 count++;
312 }
313 }
314 }
315 }
316 }
317 }
318
319 // printf("prova2\n");
320 if(ii==0 && COBT[0]<x[0]){
321 // printf("\n%f-(%f-%lld) div 10 \n",x[0],x[maxev-minev-1],COBT[0]);
322 xMIN=x[0]-(x[maxev-minev-1]-COBT[0])/10;
323 xMAX=x[maxev-minev-1]+(x[maxev-minev-1]-COBT[0])/10;
324 // printf("\nxMIN=%f\txMAX=%f\n",xMIN,xMAX);
325 }
326 else{
327 // printf("\n%f\t%f \n",x[0],x[maxev-minev-1]);
328 xMIN=x[0]-(x[maxev-minev-1]-x[0])/10;
329 xMAX=x[maxev-minev-1]+(x[maxev-minev-1]-x[0])/10;
330 // printf("\nxMIN=%f\txMAX=%f\n",xMIN,xMAX);
331 if(xMIN<0) xMIN=0;
332 }
333
334
335 for (Int_t i=0; i<12 ; i++){
336
337 Double_t y[maxev-minev],yb[maxev-minev];
338 for(Int_t v=0;v<maxev-minev;v++){
339 y[v]=yyd[v][i];
340 yb[v]=yyb[v][i];
341 }
342
343 // printf("prova3\n");
344 if((maxev-minev)>1000){
345 perc=(count*100)/(maxev-minev);
346 if(perc>10) pad[i][ii]->SetFillColor(2);
347 else pad[i][ii]->SetFillColor(10);
348 }
349 else{
350 if(count>=100) pad[i][ii]->SetFillColor(2);
351 else pad[i][ii]->SetFillColor(10);
352 }
353
354 oss<<"DSP "<<i+1;
355 DataTimeCanv[ii]->cd();
356 pad[i][ii]->SetFrameFillColor(10);
357 pad[i][ii]->Draw();
358 pad[i][ii]->cd();
359 dataletime[i][ii]= new TGraph((maxev-minev),x,y);
360 // printf("\nentries nel graph= %d\n",maxev-minev);
361 dataletime[i][ii]->SetTitle(oss.str().c_str());
362 dataletime[i][ii]->GetXaxis()->SetTitle("OBT (ms)");
363 dataletime[i][ii]->GetXaxis()->CenterTitle();
364 // dataletime[i][ii]->GetXaxis()->SetRangeUser(min,max);
365 dataletime[i][ii]->GetXaxis()->SetRangeUser(xMIN,xMAX);
366 dataletime[i][ii]->GetYaxis()->SetTitle("datalength (Word 13 bit)");
367 dataletime[i][ii]->GetYaxis()->CenterTitle();
368 if(i==6) dataletime[i][ii]->GetYaxis()->SetRangeUser(0,4500);
369 else dataletime[i][ii]->GetYaxis()->SetRangeUser(0,3500);
370 dataletime[i][ii]->SetMarkerStyle(21);
371 if((maxev-minev)<50) dataletime[i][ii]->SetMarkerSize(0.5);
372 else dataletime[i][ii]->SetMarkerSize(0.3);
373 dataletime[i][ii]->SetMarkerColor(4);
374 dataletime[i][ii]->Draw("ap");
375
376 // printf("prova4\n");
377
378 if((maxev-minev)>1000 && perc>10){
379 dataletime1[i][ii]= new TGraph(count,xb,yb);
380 dataletime1[i][ii]->SetMarkerStyle(21);
381 if((maxev-minev)<50) dataletime1[i][ii]->SetMarkerSize(0.5);
382 else dataletime1[i][ii]->SetMarkerSize(0.3);
383 dataletime1[i][ii]->SetMarkerColor(2);
384 dataletime1[i][ii]->Draw("psame");
385 }
386 else if((maxev-minev)<1000 && count>=100){
387 dataletime1[i][ii]= new TGraph(count,xb,yb);
388 dataletime1[i][ii]->SetMarkerStyle(21);
389 if((maxev-minev)<50) dataletime1[i][ii]->SetMarkerSize(0.5);
390 else dataletime1[i][ii]->SetMarkerSize(0.3);
391 dataletime1[i][ii]->SetMarkerColor(2);
392 dataletime1[i][ii]->Draw("psame");
393 }
394 li.SetLineColor(1);
395 li.SetLineStyle(1);
396 li.SetLineWidth(1);
397 // if(i!=6) li.DrawLine(min,750,max,750);
398 // else li.DrawLine(min,1500,max,1500);
399 // li.DrawLine(min,3075,max,3075);
400 if(i!=6) li.DrawLine(xMIN,750,xMAX,750);
401 else li.DrawLine(xMIN,1500,xMAX,1500);
402 li.DrawLine(xMIN,3075,xMAX,3075);
403
404 li.SetLineColor(12);
405 li.SetLineStyle(4);
406 li.SetLineWidth(1);
407 for(Int_t j=countHOBT;j<neventH;j++){
408 if(HOBT[j]<HOBT[j-1] && j!=countHOBT || HOBT[j]>x[maxev-minev-1]){
409 if(i==11) countHOBT=j;
410 break;
411 }
412 else{
413 if( HOBT[j]>xMIN && HOBT[j]<x[maxev-minev-1]){
414 if(i==6) li.DrawLine(HOBT[j],0.,HOBT[j],4500.);
415 else li.DrawLine(HOBT[j],0.,HOBT[j],3500.);
416 if(trk_cal_us[j]==104){
417 calus<<"D";
418 t2->SetTextColor(6);
419 if(i==6) t2->DrawLatex(HOBT[j],4350.,calus.str().c_str());
420 else t2->DrawLatex(HOBT[j],3350.,calus.str().c_str());
421 calus.str("");
422 }
423 else{
424 calus<<"O";
425 t2->SetTextColor(3);
426 if(i==6) t2->DrawLatex(HOBT[j],4350.,calus.str().c_str());
427 else t2->DrawLatex(HOBT[j],3350.,calus.str().c_str());
428 calus.str("");
429 }
430 }
431 }
432 }
433 // if(ii==0){
434 // if(i==6) ar.DrawArrow(COBT[0],1700.,COBT[0],2700.,0.01,"<");
435 // else ar.DrawArrow(COBT[0],1000.,COBT[0],2000.,0.01,"<");
436 // }
437 for(Int_t j=countCOBT;j<neventC;j++){
438 if(COBT[j]<COBT[j-1] && j!=countCOBT || COBT[j]>x[maxev-minev-1]){
439 if(i==11) countCOBT=j;
440 break;
441 }
442 else{
443 if( COBT[j]>xMIN && COBT[j]<x[maxev-minev-1]){
444 if(i==6) ar.DrawArrow(COBT[j],1700.,COBT[j],2700.,0.01,"<");
445 else ar.DrawArrow(COBT[j],1000.,COBT[j],2000.,0.01,"<");
446 }
447 }
448 }
449
450 oss.str("");
451 DataTimeCanv[ii]->Update();
452 }
453
454 minev=maxev;
455 printf("\ncountrun=%d\n",ii);
456 if(maxev==maxevent-1) {
457 countrun=ii+1;
458 break;
459 }
460 }
461 printf("... end of packets. \n");
462 //*************************************************************************
463 // Save output Files
464 //*************************************************************************
465 stringstream nom1,nom2,nom3;
466
467 for(Int_t fl=0;fl<countrun;fl++){
468 if(countrun==1){
469 nom1<<ffile<<"_FTrkQLook_BASIC."<<outfile.Data();
470 DataTimeCanv[fl]->Print(out+nom1.str().c_str());
471 nom1.str("");
472 }
473
474 if(countrun>=2){
475 if(!strcmp(outfile.Data(),"ps") || !strcmp(outfile.Data(),"pdf")){
476 nom1.str("");
477 nom2.str("");
478 nom3.str("");
479 nom1<<ffile<<"_FTrkQLook_BASIC.ps(";
480 nom2<<ffile<<"_FTrkQLook_BASIC.ps";
481 nom3<<ffile<<"_FTrkQLook_BASIC.ps)";
482 if(fl==0) DataTimeCanv[fl]->Print(out+nom1.str().c_str(),"portrait");
483 else if(fl==countrun-1) DataTimeCanv[fl]->Print(out+nom3.str().c_str(),"portrait");
484 else DataTimeCanv[fl]->Print(out+nom2.str().c_str(),"portrait");
485
486 }
487 else{
488 nom1.str("");
489 nom1<<ffile<<"_FTrkQLook_BASIC-pag"<<fl+1<<"."<<outfile.Data();
490 DataTimeCanv[fl]->Print(out+nom1.str().c_str());
491 }
492 }
493 }
494
495 if(!strcmp(outfile.Data(),"pdf") && countrun>=2){
496 stringstream com;
497 com<<"ps2pdf13 "<<out<<ffile<<"_FTrkQLook_BASIC.ps "<<out<<ffile<<"_FTrkQLook_BASIC.pdf";
498 system(com.str().c_str());
499 printf("\n---> ps file converted in pdf format!\n");
500 com.str("");
501 com<<"rm -f "<<out<<ffile<<"_FTrkQLook_BASIC.ps";
502 system(com.str().c_str());
503 printf("---> ps file removed!\n\n");
504 com.str("");
505 }
506
507 datafile->Close();
508 gROOT->Reset();
509 return;
510 }

  ViewVC Help
Powered by ViewVC 1.1.23