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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.10 - (hide annotations) (download)
Thu Jul 13 10:13:37 2006 UTC (18 years, 5 months ago) by pam-fi
Branch: MAIN
CVS Tags: v1r06
Changes since 1.9: +1 -1 lines
release v1r06

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

  ViewVC Help
Powered by ViewVC 1.1.23