/[PAMELA software]/quicklook/QLflightTmtc_Header/HeaderScan.cpp
ViewVC logotype

Annotation of /quicklook/QLflightTmtc_Header/HeaderScan.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.1 - (hide annotations) (download)
Fri Jun 16 10:11:54 2006 UTC (18 years, 5 months ago) by pam-rm2
Branch: MAIN
Branch point for: QLflightTmtc_Header
Initial revision

1 pam-rm2 1.1 /**
2     * Header Scan
3     * Author Nagni
4     * version 1.0
5     *
6     * Version 1.1 - 28 December 2004
7     * If outList does not exist the function exit to prompt
8     * If nevents = 0 the function exit to promt
9     *
10     * Version 1.2 - 3 January 2005
11     * Two canvases are created to see the graphs better
12     *
13     * Version 1.25 - 13 January 2005
14     * Changed Int_t to Float because of variable range size
15     * (UInt_t has been excluded beacuse of uncompatibility with TGraph)
16     *
17     * version 1.3 - 22 February 2005 - Nagni
18     * For compatibility with batch mode excution:
19     * 1) Added "include <iostream>" and "using namespace std"
20     * 2) Removed gROOT->Reset()
21     *
22     * Version 1.4
23     * Date 08 March 2005
24     * Added "format" parameter for saving the produced image in various formats
25     * (for a complete list of types refer to TPad::SaveAs method)
26     *
27     * Version 1.5
28     * Date 09 February 2006 - Marcelli
29     * Update to work with new Yoda output
30     *
31     * Version 1.6
32     * Date 27 February 2006 - Marcelli
33     * For compilation:
34     * Added function "int main(int argc, char* argv[])"
35     *
36     *
37     * Description: This script creates two canvases with five pads. The first pad shows packetID variable (for all packets) vs. OBT.
38     * The second pad shows the number of physic packets vs. OBT. The third pad shows the lenght of Physic packets (byte) vs. OBT.
39     * The fourth pad shows the packetcounter of physic packets vs. OBT. The fifth pad shows PacketCounter vs. File Offset.
40     *
41     * Parameters:
42     * TSTring base - the path to the root directory for the specific Pamela unpack session
43     * There is no default value, without this input the program will not run
44     * TString outDir - the path where to save the output image (Default = ./)
45     * TString format - the format which will be used fo rsave the produced images (Default = "jpg")
46     */
47    
48    
49    
50     #include <sstream>
51     #include <iostream>
52     #include "TString.h"
53     #include "TStyle.h"
54     #include "TFile.h"
55     #include "TList.h"
56     #include "TTree.h"
57     #include "TObjString.h"
58     #include "TCanvas.h"
59     #include "TGraph.h"
60     #include "TH1F.h"
61     #include "TGaxis.h"
62     #include "EventHeader.h"
63     #include "PscuHeader.h"
64    
65    
66     using namespace std;
67    
68     void HeaderScan(TString base, TString outDir, TString format)
69     {
70    
71     TList *list = new TList;
72     Int_t numkey;
73     TObject *key = new TObject;
74     const char *name;
75     TTree* tr = new TTree;
76     Long64_t totevents=0, totphysevents=0;
77     Float_t id;
78     Long64_t nevents=0;
79    
80     //------- load root file --------------
81    
82     TFile *file = new TFile(base.Data());
83    
84     if (!file){
85     printf("No such file in the directory has been found");
86     return;
87     }
88     if (outDir == "" ) outDir = ".";
89    
90     list = file->GetListOfKeys();
91    
92     numkey = file->GetNkeys();
93    
94     pamela::EventHeader *eh=0;
95     pamela::PscuHeader *ph=0;
96    
97     ///-----to know the total number f events end of physics events----//////
98    
99     for (Int_t i=0; i<numkey; i++){
100     key = list->At(i);
101     name=(key->GetName());
102     //cout<<name<<"\n";
103     tr = (TTree*)file->Get(name);
104     if (tr->IsZombie()) continue;
105     nevents = tr->GetEntries();
106     totevents+=nevents;
107     tr->SetBranchAddress("Header", &eh);
108     for (Int_t j = 0; j < nevents; j++){
109     tr->GetEntry(j);
110     ph = eh->GetPscuHeader();
111     if(ph->GetPacketId1() == 0x10) ++totphysevents;
112     }
113    
114     }
115    
116     const Long64_t totalevents=totevents;
117     const Long64_t totalphysevents=totphysevents;
118    
119     /////////////////////////////////////////////////
120    
121     TString filename = ((TObjString*)base.Tokenize('/')->Last())->GetString();
122     filename = ((TObjString*)filename.Tokenize('.')->First())->GetString();
123     stringstream oss;
124     oss.str("");
125     oss << filename.Data();
126    
127     Float_t obt[totalevents], pcktId[totalevents], pcktLenght[totalphysevents], pcktCounter[totalphysevents], offset[totalevents], allCounter[totalevents];
128     Float_t obtphysevents=0;
129    
130    
131     totphysevents=0;
132     totevents=0;
133     for (Int_t i=0; i<numkey; i++){
134     key = list->At(i);
135     name=(char *)(key->GetName());
136     tr = (TTree*)file->Get(name);
137     if (tr->IsZombie()) continue;
138     tr->SetBranchAddress("Header", &eh);
139     nevents = tr->GetEntries();
140     for (Int_t j = 0; j < nevents; j++){
141     tr->GetEntry(j);
142     ph = eh->GetPscuHeader();
143     obt[j+totevents]=ph->GetOrbitalTime() ;
144     pcktId[j+totevents]=ph->GetPacketId1();
145     offset[j+totevents]=ph->GetFileOffset();
146     if(pcktId[j+totevents]==16){
147     pcktLenght[totphysevents]=ph->GetPacketLenght();
148     pcktCounter[totphysevents]=ph->GetCounter();
149     totphysevents=totphysevents+1;
150     }
151     offset[j+totevents]=ph->GetFileOffset();
152     allCounter[j+totevents]= ph->GetCounter();
153     }
154     totevents=totevents+nevents;
155     }
156    
157    
158     Float_t mintime=obt[0], maxtime=obt[0], minlen=pcktLenght[0], maxlen=pcktLenght[0], mincount=pcktCounter[0], maxcount=pcktCounter[0];
159    
160    
161     for (Int_t t=0; t<totalevents; t++){
162     if(obt[t]<mintime) mintime=obt[t];
163     if(obt[t]>maxtime) maxtime=obt[t];
164     }
165     for (Int_t t=0; t<totalphysevents; t++){
166     if(pcktLenght[t]<=minlen) minlen=pcktLenght[t];
167     if(pcktLenght[t]>=maxlen) maxlen=pcktLenght[t];
168     if(pcktCounter[t]<=mincount) mincount=pcktCounter[t];
169     if(pcktCounter[t]>=maxcount) maxcount=pcktCounter[t];
170     }
171    
172    
173     Float_t nbintime = (maxtime-mintime)/1000;
174     Float_t nbinlen = (maxlen-minlen)/100;
175     Float_t nbincount = (maxcount-mincount)/100;
176    
177    
178     ///---------------------------TO GRAPH---------------------------------------------///
179     TCanvas *finalCanv1 = new TCanvas("Header_1", base, 1280, 1024);
180     TCanvas *finalCanv2 = new TCanvas("Header_2", base, 1280, 1024);
181     finalCanv1->Divide(3);
182     finalCanv1->SetFillColor(10);
183     finalCanv2->Divide(2);
184     finalCanv2->SetFillColor(10);
185    
186    
187     TGraph *grPcktId1 = new TGraph(totalevents, obt, pcktId);
188     TGraph *grPcktLenght = new TGraph(totalevents, obt, pcktLenght);
189     TGraph *grPcktCounter = new TGraph(totalevents, obt, pcktCounter);
190     TGraph *grFileOffset = new TGraph(totalevents, offset, allCounter);
191    
192     TGaxis::SetMaxDigits(3);
193    
194     //-----canvas 1-------------------------------///
195     finalCanv1->cd(1);
196    
197     stringstream oss1;
198     oss1.str("");
199     oss1 << "OBT (ms) with t0 = " << mintime << "ms";
200    
201     gStyle->SetTitleH(0.06);
202     grPcktId1->SetTitle(oss.str().c_str());
203     grPcktId1->GetXaxis()->SetTitle(oss1.str().c_str());
204     grPcktId1->GetXaxis()->CenterTitle();
205     grPcktId1->GetXaxis()->SetLabelSize(0.03);
206     grPcktId1->GetYaxis()->SetTitle("Packet ID");
207     grPcktId1->GetYaxis()->CenterTitle();
208     grPcktId1->SetMarkerSize(4);
209     grPcktId1->Draw("AP");
210    
211    
212     finalCanv1->cd(2);
213     oss1.str("");
214     oss1 << "OBT (min) with t0 = " << mintime << "ms";
215    
216     TH1F *h1 = new TH1F ("h1", oss.str().c_str(), (int)(nbintime/60), mintime, maxtime);
217    
218     for (Int_t i=0; i<numkey; i++){
219     key = list->At(i);
220     name=(char *)(key->GetName());
221     tr = (TTree*)file->Get(name);
222     if (tr->IsZombie()) continue;
223     tr->SetBranchAddress("Header", &eh);
224     nevents = tr->GetEntries();
225     for (Int_t j = 0; j < nevents; j++){
226     tr->GetEntry(j);
227     ph = eh->GetPscuHeader();
228     if((ph->GetPacketId1()) == 16){
229     obtphysevents=ph->GetOrbitalTime();
230     h1->Fill(obtphysevents);
231     }
232     }
233     }
234    
235     h1->SetMarkerStyle(6);
236     h1->GetXaxis()->SetTitle(oss1.str().c_str());
237     h1->GetXaxis()->CenterTitle();
238     h1->GetXaxis()->SetLabelSize(0.03);
239     h1->GetYaxis()->SetTitle("number of Physic packets");
240     h1->GetYaxis()->CenterTitle();
241     h1->Draw();
242    
243    
244     finalCanv1->cd(3);
245     oss1.str("");
246     oss1 << "OBT (ms) with t0 = " << mintime << "ms";
247    
248     grPcktLenght->SetTitle(oss.str().c_str());
249     grPcktLenght->GetXaxis()->SetTitle(oss1.str().c_str());
250     grPcktLenght->GetXaxis()->CenterTitle();
251     grPcktLenght->GetXaxis()->SetLabelSize(0.03);
252     grPcktLenght->GetYaxis()->SetTitle("Lenght of Physic packets (byte)");
253     grPcktLenght->GetYaxis()->CenterTitle();
254     grPcktLenght->GetYaxis()->SetLabelSize(0.03);
255     grPcktLenght->SetMarkerSize(4);
256     grPcktLenght->Draw("AP");
257    
258    
259     finalCanv1->cd(2);
260     h1->Draw();
261    
262     finalCanv1->Update();
263    
264    
265     ///---------canvas 2-----------------------------//
266     finalCanv2->cd(1);
267    
268     grPcktCounter->SetTitle(oss.str().c_str());
269     grPcktCounter->GetXaxis()->SetTitle(oss1.str().c_str());
270     grPcktCounter->GetXaxis()->SetTitleSize(0.04);
271     grPcktCounter->GetXaxis()->CenterTitle();
272     grPcktCounter->GetXaxis()->SetLabelSize(0.03);
273     grPcktCounter->GetYaxis()->SetTitle("PacketCounter of Physic packets");
274     grPcktCounter->GetYaxis()->SetTitleSize(0.04);
275     grPcktCounter->GetYaxis()->CenterTitle();
276     grPcktCounter->GetYaxis()->SetLabelSize(0.03);
277     grPcktCounter->SetMarkerSize(4);
278     grPcktCounter->SetMinimum(mincount);
279     grPcktCounter->SetMaximum(maxcount);
280     grPcktCounter->Draw("AP");
281    
282     finalCanv2->cd(2);
283    
284     grFileOffset->SetTitle(oss.str().c_str());
285     grFileOffset->GetXaxis()->SetTitle("File Offset");
286     grFileOffset->GetXaxis()->CenterTitle();
287     grFileOffset->GetXaxis()->SetTitleSize(0.04);
288     grFileOffset->GetXaxis()->SetLabelSize(0.03);
289     grFileOffset->GetYaxis()->SetTitle("Packet counter");
290     grFileOffset->GetYaxis()->CenterTitle();
291     grFileOffset->GetYaxis()->SetTitleSize(0.04);
292     grFileOffset->GetYaxis()->SetLabelSize(0.03);
293     grFileOffset->SetMarkerSize(4);
294     grFileOffset->Draw("AP");
295    
296     finalCanv1->Update();
297    
298    
299     //-------to save---------------------------///
300     oss.str("");
301     oss1.str("");
302    
303     oss << outDir.Data() << filename.Data();
304     oss1 << outDir.Data() << filename.Data();
305    
306     oss << "_HeaderScan_1." << format.Data();
307     oss1 << "_HeaderScan_2." << format.Data();
308    
309     finalCanv1->SaveAs(oss.str().c_str());
310     finalCanv2->SaveAs(oss1.str().c_str());
311    
312    
313     file->Close();
314    
315     }
316    
317    
318    
319     int main(int argc, char* argv[]){
320     TString path;
321     TString outDir = "./";
322     TString format = "jpg";
323    
324     if (argc < 2){
325     printf("You have to insert at least the file to analyze \n");
326     printf("Try '--help' for more information. \n");
327     exit(1);
328     }
329    
330     if (!strcmp(argv[1], "--help")){
331     printf( "Usage: HeaderScan FILE [OPTION] \n");
332     printf( "\t --help Print this help and exit \n");
333     printf( "\t -outDir[path] Path where to put the output [default ./] \n");
334     printf( "\t -format[jpg|ps|gif] Format for output files [default 'jpg'] \n");
335     exit(1);
336     }
337    
338    
339     path=argv[1];
340    
341     for (int i = 2; i < argc; i++){
342    
343     if (!strcmp(argv[i], "-outDir")){
344     if (++i >= argc){
345     printf( "-outDir needs arguments. \n");
346     printf( "Try '--help' for more information. \n");
347     exit(1);
348     }
349     else{
350     outDir = argv[i];
351     continue;
352     }
353     }
354    
355     if (!strcmp(argv[i], "-format")){
356     if (++i >= argc){
357     printf( "-format needs arguments. \n");
358     printf( "Try '--help' for more information. \n");
359     exit(1);
360     }
361     else{
362     format = argv[i];
363     continue;
364     }
365     }
366    
367    
368     }
369    
370     HeaderScan(argv[1], outDir, format);
371    
372     }

  ViewVC Help
Powered by ViewVC 1.1.23