/[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.5 - (hide annotations) (download)
Mon Jun 5 14:23:29 2006 UTC (18 years, 6 months ago) by pam-fi
Branch: MAIN
Changes since 1.4: +3 -17 lines
*** empty log message ***

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

  ViewVC Help
Powered by ViewVC 1.1.23