/[PAMELA software]/quicklook/tracker/flight/src/TrkQLook_BASIC.cpp
ViewVC logotype

Annotation of /quicklook/tracker/flight/src/TrkQLook_BASIC.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.1 - (hide annotations) (download)
Tue Mar 14 16:09:30 2006 UTC (18 years, 10 months ago) by pam-fi
Branch: MAIN
Branch point for: trk-QLook
Initial revision

1 pam-fi 1.1 #include <utils/yodaUtility.h>
2     #include "TrkFunctions.cpp"
3     #include <iostream>
4     #include <TPaveText.h>
5     #include <TLatex.h>
6     #include <TCanvas.h>
7     #include <TGraph.h>
8     #include <PscuHeader.h>
9     #include <EventHeader.h>
10     #include <RunHeaderEvent.h>
11     #include <TStyle.h>
12     #include <TObjString.h>
13     #include <physics/tracker/TrackerEvent.h>
14    
15     /**
16     * TrkQLook_BASIC.cpp
17     *
18     * autor: D.Fedele
19     * version 1.0
20     * Parameters:
21     * base - the path to the root directory for the specific
22     * Pamela unpack session
23     * fromevent - first event to analyse
24     * toevent - last event to analyse
25     * outdir - total path of output file
26     * outfile - extension of output file (pdf,ps,gif,jpg)
27     *
28     */
29    
30     void TrkQLook_BASIC(TString file,Int_t fromevent,Int_t toevent, TString outdir,TString outfile)
31     {
32     //
33     // obtain information about the data file and select the output dir
34     const string filepath=file.Data();
35     Int_t dwpos = filepath.find("DW_");
36     Int_t dwpos1 = filepath.find(".root");
37     TString fpath=(filepath.c_str());
38     TString base,ffile ;
39     stringcopy(ffile,fpath,dwpos,dwpos1);
40     stringcopy(base,fpath,0,dwpos);
41    
42     TString out;
43     if(outdir.Length()==0){
44     out = base;
45     }else{
46     out = outdir;
47     }
48    
49     //
50     // inizialise the variables and open the file
51     pamela::tracker::TrackerEvent *te=0;
52     pamela::EventHeader *eh=0,*eH=0;
53     pamela::RunHeaderEvent *reh=0;
54     pamela::PscuHeader *ph=0,*pH=0;
55    
56     TFile *datafile = new TFile(file);
57     TTree *otr = (TTree*)datafile->Get("RunHeader");
58     otr->SetBranchAddress("Header",&eH);
59     otr->SetBranchAddress("RunHeader",&reh);
60     TTree *tr = (TTree*)datafile->Get("Physics");
61     tr->SetBranchAddress("Tracker",&te);
62     tr->SetBranchAddress("Header",&eh);
63    
64    
65     Long64_t nevent = tr->GetEntries();
66     Long64_t neventH = otr->GetEntries();
67     Int_t minevent=0;
68     Int_t maxevent=0;
69     // cout<<"Total n.events "<<nevent<<" e neventH= "<<neventH<<endl;
70    
71     if (nevent<=0){
72     datafile->Close();
73     return;
74     }
75     if ( fromevent > toevent && toevent>0 ){
76     printf("It must be fromevent < toevent \n");
77     return;
78     }
79     if ( fromevent > nevent || fromevent < 0 ) {
80     printf("You can choose fromevent between 0 (all) and %i \n",nevent);
81     return;
82     }
83     if ( toevent > nevent || toevent < 0 ) {
84     printf("You can choose toevent between 0 (all) and %i \n",nevent);
85     return;
86     }
87     if ( fromevent == 0 ) {
88     minevent = 0;
89     maxevent = nevent;
90     } else {
91     minevent = fromevent;
92     if ( toevent > 0 ){
93     maxevent = toevent+1;
94     } else if (toevent > nevent) {
95     maxevent = nevent;
96     } else {
97     maxevent = toevent+1;
98     }
99     nevent=maxevent-minevent ;
100     }
101    
102     //
103     // other variables definitions
104     stringstream oss,fromfile,isfile;
105     const Int_t size=nevent;
106     Int_t dsp=0,count=0,trk_cal_us[50];
107     Double_t perc;
108     Double_t yd[size*12], x[size+1];
109     Double_t yyd[size+1],yyb[size], xb[size];
110     ULong64_t HOBT[50];
111     TGraph *dataletime,*dataletime1;
112     TPaveText *pt1;
113     TPad *pt;
114    
115     for (Int_t i=0; i<50;i++){
116     HOBT[i]=0;
117     trk_cal_us[i]=0;
118     }
119    
120     //***************************************************************************************
121     // LOOP on each event
122     //***************************************************************************************
123    
124     //
125     // information about trk_calib_used
126     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     }
132    
133    
134     //
135     // information about datalength
136     x[0]=(Double_t)HOBT[0];
137     for (Int_t ev=minevent; ev<maxevent; ev++){
138     tr->GetEntry(ev);
139     ph = eh->GetPscuHeader();
140    
141     x[(ev-minevent)+1]= ph->GetOrbitalTime();
142     dsp=0;
143     for(Int_t i=0; i<12; i++){
144     dsp=te->DSPnumber[i]-1;
145     yd[(ev-minevent)*12+dsp]= te->DATAlength[i];
146     }
147     };
148     datafile->Close();
149    
150    
151     //****************************************************************************************
152     //Output figure
153     //****************************************************************************************
154    
155     TCanvas *DataTimeCanv=new TCanvas("Tracker_Detector_Report_","",900,1200);
156     DataTimeCanv->SetFillColor(10);
157     DataTimeCanv->Range(0,0,100,100);
158     fromfile<<"TrackerQLookBasic File: "<<ffile;
159     isfile<<"DATALENGTH vs. ORBITALTIME ";
160    
161     TLatex *t=new TLatex();
162     t->SetTextFont(32);
163     t->SetTextColor(1);
164     t->SetTextAlign(12);
165     t->SetTextSize(0.02);
166     t->DrawLatex(2.,98.7,fromfile.str().c_str());
167     TLatex *t1=new TLatex();
168     t1->SetTextFont(32);
169     t1->SetTextColor(1);
170     t1->SetTextAlign(12);
171     t1->SetTextSize(0.02);
172     t1->DrawLatex(70.,98.7,isfile.str().c_str());
173     fromfile.str("");
174     isfile.str("");
175    
176    
177     //*************************************************************************************
178     //book pads and histos
179     //***************************************************************************************
180     gStyle->SetLabelSize(0.06,"x");
181     gStyle->SetLabelSize(0.06,"y");
182     gStyle->SetStatFontSize(0.075);
183     gStyle->SetOptStat(1110);
184     gStyle->SetFillColor(10);
185     gStyle->SetTitleFontSize(0.1);
186     gStyle->SetTitleOffset(0.8,"y");
187     gStyle->SetTitleOffset(0.9,"x");
188     gStyle->SetTitleSize(0.06,"y");
189     gStyle->SetTitleSize(0.055,"x");
190    
191     TPad *pad[12]; //pad for histos
192     Double_t posy = 0.95; // up y-coord - top pads
193     Double_t hpad = 0; // pad height
194     Double_t posx1=0; // left x-coord - pad column
195     Double_t posx0=0; // x-coord - column division
196     Double_t wrel = 0; // relative x size of first sub-column
197     Double_t marg = 0.004; // margin among pads
198    
199     hpad = (posy-marg*11)/6;
200     wrel = (1-marg*4)/2;
201     stringstream title;
202     stringstream hid;
203    
204     for(Int_t n = 0; n<12; n++) {
205     if ( (n+1)%2==1 ) {
206     if(n>1) posy = posy-(marg*2+hpad);
207     posx1 = marg;
208     posx0 = posx1 + wrel;
209     }
210     else {
211     posx1 = posx0 + 2*marg;
212     posx0 = posx1 + wrel;
213     }
214    
215     /* -----------> pad for histograms */
216     pad[n] = new TPad("pad"," ",posx1,posy-hpad,posx0,posy,18,0,0);
217     };
218    
219     TLine li;
220     li.SetLineColor(1);
221     li.SetLineStyle(1);
222     li.SetLineWidth(1);
223     //**********************************************************************************
224     // Fill Graphs and Histos
225     //**********************************************************************************
226     stringstream calus;
227     TLatex *t2=new TLatex();
228     t2->SetTextFont(32);
229     t2->SetTextColor(1);
230     t2->SetTextAlign(13);
231     t2->SetTextSize(0.08);
232    
233     for (Int_t i=0; i<12 ; i++){
234     perc=0;
235     count=0;
236     yyd[0]=-10000.;
237     for (Int_t ev=minevent; ev<maxevent; ev++){
238     yyd[(ev-minevent)+1]=yd[12*(ev-minevent)+i];
239     if(i==6){
240     if(yyd[(ev-minevent)+1]>1500){
241     if(yyd[(ev-minevent)+1]<3075){
242     yyb[count]= yyd[(ev-minevent)+1];
243     xb[count]= x[(ev-minevent)+1];
244     count++;
245     }
246     }
247     }
248     else{
249     if(yyd[(ev-minevent)+1]>750){
250     if(yyd[(ev-minevent)+1]<3075){
251     yyb[count]= yyd[(ev-minevent)+1];
252     xb[count]= x[(ev-minevent)+1];
253     count++;
254     }
255     }
256     }
257     }
258    
259     if((maxevent-minevent)>100) perc=(count*100)/(maxevent-minevent);
260     else perc=(count*10)/(maxevent-minevent);
261    
262     Double_t min,max;
263     if((maxevent-minevent)<100){
264     min=x[0];
265     max=x[size]*1.0001;
266     }
267     else{
268     min=x[0];
269     max=x[size-1]*1.0005;
270     }
271     oss<<"DSP "<<i+1;
272     DataTimeCanv->cd();
273     if(perc>1) pad[i]->SetFillColor(2);
274     else pad[i]->SetFillColor(10);
275     pad[i]->SetFrameFillColor(10);
276     pad[i]->Draw();
277     pad[i]->cd();
278     dataletime= new TGraph((maxevent-minevent)+1,x,yyd);
279     dataletime->SetTitle(oss.str().c_str());
280     dataletime->GetXaxis()->SetTitle("OBT (ms)");
281     dataletime->GetXaxis()->CenterTitle();
282     dataletime->GetXaxis()->SetRangeUser(min,max);
283     dataletime->GetYaxis()->SetTitle("datalength (Word 13 bit)");
284     dataletime->GetYaxis()->CenterTitle();
285     dataletime->GetYaxis()->SetRangeUser(0,3500);
286     dataletime->SetMarkerStyle(21);
287     if((maxevent-minevent)<50) dataletime->SetMarkerSize(0.5);
288     else dataletime->SetMarkerSize(0.3);
289     dataletime->SetMarkerColor(1);
290     dataletime->Draw("ap");
291     if(perc>1){
292     dataletime1= new TGraph(count,xb,yyb);
293     dataletime1->SetMarkerStyle(21);
294     if((maxevent-minevent)<50) dataletime1->SetMarkerSize(0.5);
295     else dataletime1->SetMarkerSize(0.3);
296     dataletime1->SetMarkerColor(2);
297     dataletime1->Draw("psame");
298     }
299     li.SetLineColor(1);
300     li.SetLineStyle(1);
301     li.SetLineWidth(1);
302     if(i!=6) li.DrawLine(min,750,max,750);
303     else li.DrawLine(min,1500,max,1500);
304     li.DrawLine(min,3075,max,3075);
305     for(Int_t j=0;j<neventH;j++){
306     li.SetLineColor(12);
307     li.SetLineStyle(4);
308     li.SetLineWidth(1);
309     li.DrawLine(HOBT[j],0.,HOBT[j],3500.);
310     if(trk_cal_us[j]==104){
311     calus<<"D";
312     t2->SetTextColor(6);
313     t2->DrawLatex(HOBT[j],3350.,calus.str().c_str());
314     calus.str("");
315     }
316     else{
317     calus<<trk_cal_us[j];
318     t2->SetTextColor(3);
319     t2->DrawLatex(HOBT[j],3350.,calus.str().c_str());
320     calus.str("");
321     }
322     }
323     oss.str("");
324     };
325    
326     DataTimeCanv->Update();
327    
328     //*************************************************************************
329     // Save output Files
330     //*************************************************************************
331     stringstream out1;
332    
333     // if(!strcmp(outfile.Data(),"pdf")){
334     // stringstream com;
335     // out1<<ffile<<"_TrkQLook_BASIC.ps";
336     // DataTimeCanv->Print(out+out1.str().c_str());
337     // com<<"ps2pdf13 "<<out<<ffile<<"_TrkQLook_BASIC.ps "<<out<<ffile<<"_TrkQLook_BASIC.pdf";
338     // system(com.str().c_str());
339     // com.str("");
340     // com<<"rm -f "<<out<<ffile<<"_TrkQLook_BASIC.ps ";
341     // system(com.str().c_str());
342     // com.str("");
343     // }
344     // else{
345     out1<<ffile<<"_TrkQLook_BASIC."<<outfile.Data();
346     DataTimeCanv->Print(out+out1.str().c_str());
347     out1.str("");
348     // }
349    
350     dataletime->Delete();
351     delete DataTimeCanv;
352    
353     gROOT->Reset();
354     // printf("... end of packets. \n");
355     return;
356     }

  ViewVC Help
Powered by ViewVC 1.1.23