/[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.8 - (hide annotations) (download)
Thu Aug 31 14:28:32 2006 UTC (18 years, 3 months ago) by pam-rm2
Branch: MAIN
CVS Tags: v2r04
Changes since 1.7: +229 -135 lines
modificato HeaderScan: unica canvas, identificazione run

1 pam-rm2 1.1 /**
2 pam-rm2 1.4 * Header_Scan
3 pam-rm2 1.1 * 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 pam-rm2 1.2 *
47     *
48     * Version 1.7
49     * Date 16 June 2006 - Malvezzi
50     *
51     * Description of changes:
52     * Implementation of the case: numebr of events <= 0.
53     * Remove graph "grPcktId1"; see PacketScan for the same information.
54     * Fixed bugs: for a large namber of events is not possible to have vectors, so I have subsituted graphs with histograms
55     * or divided the graphs in more than one canvas.
56     *
57 pam-rm2 1.6 * Version 1.8
58     * Date 8 August 2006 - Malvezzi
59     *
60     * Description: changed the scale in the second and third graph of the first canvas; added a pad of text in the second canvas
61     *
62 pam-rm2 1.1 */
63    
64    
65 pam-rm2 1.2 #include <fstream>
66     #include <math.h>
67     #include "TLatex.h"
68     #include "TF1.h"
69     #include "TPaveText.h"
70     #include "TMultiGraph.h"
71 pam-rm2 1.1 #include <sstream>
72     #include <iostream>
73     #include "TString.h"
74     #include "TStyle.h"
75     #include "TFile.h"
76     #include "TList.h"
77     #include "TTree.h"
78     #include "TObjString.h"
79     #include "TCanvas.h"
80     #include "TGraph.h"
81     #include "TH1F.h"
82     #include "TGaxis.h"
83     #include "EventHeader.h"
84     #include "PscuHeader.h"
85 pam-rm2 1.8 #include "RunHeaderEvent.h"
86     #include "TPaveText.h"
87 pam-rm2 1.1
88     using namespace std;
89    
90     void HeaderScan(TString base, TString outDir, TString format)
91     {
92    
93 pam-rm2 1.8
94 pam-rm2 1.2 //------------------- Variables initilization -------------------------//
95 pam-rm2 1.8 Long64_t nevents=0, runnevents;
96 pam-rm2 1.7 ULong_t lastime, firstime, primotempo, ultimotempo, primoffset=500000000, ultimoffset;
97 pam-rm2 1.8 double obmin=0., time=0.;
98 pam-rm2 1.2 double obmax=0.;
99 pam-rm2 1.6 stringstream oss, oss1, oss2, oss3, noentries, stringa;
100 pam-rm2 1.1 //------- load root file --------------
101     TFile *file = new TFile(base.Data());
102     if (!file){
103     printf("No such file in the directory has been found");
104     return;
105     }
106    
107 pam-rm2 1.2 TTree *PhysicsTr = (TTree*)file->Get("Physics");
108     TBranch *headBr = PhysicsTr->GetBranch("Header");
109    
110     pamela::EventHeader *eh = 0;
111     pamela::PscuHeader *ph = 0;
112    
113     PhysicsTr->SetBranchAddress("Header", &eh);
114     nevents = headBr->GetEntries();
115 pam-rm2 1.8
116     TTree *RunHeadTr = (TTree*)file->Get("RunHeader"); ///run header tree
117     pamela::EventHeader *eH=0;
118     pamela::RunHeaderEvent *reh=0;
119    
120     RunHeadTr->SetBranchAddress("Header",&eH);
121     RunHeadTr->SetBranchAddress("RunHeader",&reh);
122     runnevents = RunHeadTr->GetEntries();
123    
124 pam-rm2 1.1
125     TString filename = ((TObjString*)base.Tokenize('/')->Last())->GetString();
126     filename = ((TObjString*)filename.Tokenize('.')->First())->GetString();
127 pam-rm2 1.2
128     //----------- If nevents < = 0 ---------------------------------/
129     if (nevents<=0) {
130     printf("nevents = %i \n", nevents);
131     printf(" \n");
132    
133     TCanvas *canv = new TCanvas("No entries", "No entries ", 400, 200);
134     canv->SetFillColor(10);
135     canv->cd();
136    
137     TLatex *l = new TLatex();
138     l->SetTextAlign(12);
139     l->SetTextSize(0.15);
140     l->SetTextColor(2);
141     noentries.str("");
142     noentries<< "HeaderScan_QL:";
143     l->DrawLatex(0.05, 0.7, noentries.str().c_str());
144     noentries.str("");
145     noentries<< "No Physics entries for this files";
146     l->DrawLatex(0.05, 0.5, noentries.str().c_str());
147    
148     if (outDir == "./") {
149     oss.str("");
150     oss << filename.Data() << "_HeaderScan_QL." << format.Data();
151     } else {
152     oss.str("");
153     oss << outDir.Data() << filename.Data() << "_HeaderScan_QL." << format.Data();
154     }
155 pam-rm2 1.1
156 pam-rm2 1.2 canv->Update();
157     canv->SaveAs(oss.str().c_str());
158 pam-rm2 1.1
159 pam-rm2 1.2 return;
160 pam-rm2 1.1 }
161 pam-rm2 1.2
162     //-------------- to know the max and min OBT ----------------------------//
163 pam-rm2 1.5 headBr->GetEntry(0);
164     ph = eh->GetPscuHeader();
165     firstime = ph->GetOrbitalTime();
166    
167 pam-rm2 1.8 int iii =0;
168     while(iii < nevents){
169     headBr->GetEntry(iii);
170     ph = eh->GetPscuHeader();
171     if((ph->GetOrbitalTime()) <= firstime) firstime=ph->GetOrbitalTime();
172     if((ph->GetOrbitalTime()) >= lastime) lastime=ph->GetOrbitalTime();
173     iii++;
174     }
175 pam-rm2 1.4
176 pam-rm2 1.2 //------------------------ First histogram -----------------------------------//
177 pam-rm2 1.6
178 pam-rm2 1.8 obmin=firstime;
179     obmax=lastime;
180 pam-rm2 1.4
181 pam-rm2 1.2 oss1.str("");
182 pam-rm2 1.6 oss1 << filename.Data() <<": Physics Packet per minute. Start time = " << obmin << ", End time = "<< obmax <<" ms";
183 pam-rm2 1.2 Int_t nbin = (lastime-firstime)/60000;
184 pam-rm2 1.8 //TH1F *h1 = new TH1F ("histo1", oss1.str().c_str(), nbin, obmin, obmax);
185     TH1F *h1 = new TH1F ("histo1", "" , nbin, obmin, obmax);
186    
187 pam-rm2 1.2 oss2.str("");
188 pam-rm2 1.4 oss2 << filename.Data() <<": Lenght of Physic packets";
189 pam-rm2 1.6 Int_t nbin2 =(lastime-firstime);
190 pam-rm2 1.8 //TH1F *packetLength = new TH1F ("packetLength", oss2.str().c_str(), nbin2, obmin, obmax);
191     TH1F *packetLength = new TH1F ("packetLength", "", nbin2, obmin, obmax);
192 pam-rm2 1.2
193     oss3.str("");
194 pam-rm2 1.8 oss3 << filename.Data() <<": Physics Counter vs. OBT";
195     // TH1F *packeCounter = new TH1F ("packeCounter", oss3.str().c_str(), nbin2, obmin, obmax);
196     TH1F *packeCounter = new TH1F ("packeCounter", "", nbin2, obmin, obmax);
197 pam-rm2 1.6 //----------------------------------------------------
198 pam-rm2 1.8 TCanvas *finalCanv = new TCanvas("Header", base, 1200, 1600);
199     finalCanv->Divide(1,5);
200     finalCanv->SetFillColor(10);
201    
202 pam-rm2 1.6 TPad *all2= new TPad ("","", 0, 0, 1, 1);
203     all2->SetFillColor(10);
204     TPad *all3= new TPad ("","", 0, 0, 1, 1);
205     all3->SetFillColor(10);
206     TPad *all4= new TPad ("","", 0, 0, 1, 1);
207     all4->SetFillColor(10);
208 pam-rm2 1.8 TPad *all= new TPad ("","", 0, 0, 1, 1);
209     all->SetFillColor(10);
210     TPad *all1= new TPad ("","", 0, 0, 1, 1);
211     all1->SetFillColor(10);
212     TPad *pad = new TPad("pad","pad", .80,.50,.90,.80);
213     pad->SetFillColor(10);
214 pam-rm2 1.6 //----------------------------------------------------
215 pam-rm2 1.4
216 pam-rm2 1.6 TMultiGraph *mg1 = new TMultiGraph();
217 pam-rm2 1.8
218 pam-rm2 1.6
219     TMultiGraph *mg2 = new TMultiGraph();
220 pam-rm2 1.4
221 pam-rm2 1.8
222 pam-rm2 1.6 //--------------------------------------------------------------------------
223 pam-rm2 1.8 for (Int_t l = 0; l < nevents; l++){
224 pam-rm2 1.6 headBr->GetEntry(l);
225 pam-rm2 1.2 ph = eh->GetPscuHeader();
226 pam-rm2 1.6 h1->Fill(ph->GetOrbitalTime());
227     packetLength->Fill(ph->GetOrbitalTime(),ph->GetPacketLenght());
228     packeCounter->Fill(ph->GetOrbitalTime(),ph->GetCounter());
229 pam-rm2 1.8 }
230 pam-rm2 1.4
231     //------------ First Canvas ---------------------//
232 pam-rm2 1.8 TLine li;
233     li.SetLineStyle(4);
234     li.SetLineWidth(2);
235    
236     finalCanv->cd(1);
237 pam-rm2 1.6 all2->Draw();
238     all2->cd();
239 pam-rm2 1.8
240 pam-rm2 1.2 h1->SetStats(kFALSE);
241     h1->GetXaxis()->SetTitle("OBT (ms)");
242 pam-rm2 1.1 h1->GetXaxis()->CenterTitle();
243 pam-rm2 1.2 h1->GetXaxis()->SetLabelSize(0.04);
244     h1->GetYaxis()->SetTitle("Number of events ");
245 pam-rm2 1.1 h1->GetYaxis()->CenterTitle();
246 pam-rm2 1.8 h1->GetYaxis()->SetLabelSize(0.06);
247 pam-rm2 1.2 h1->GetYaxis()->SetTitleSize(0.06);
248 pam-rm2 1.8 h1->GetYaxis()->SetTitleOffset(0.6);
249 pam-rm2 1.1 h1->Draw();
250 pam-rm2 1.8 for (Int_t l = 0; l < runnevents; l++){
251     RunHeadTr->GetEntry(l);
252     ph = eH->GetPscuHeader();
253     int ws= reh->RM_ACQ_SETTING_MODE;
254     int id = ph->GetPacketId1();
255     Int_t obt = ph->GetOrbitalTime();
256     if (ws==1){
257     li.SetLineColor(3);
258     li.DrawLine(obt,0,obt,h1->GetMaximum());
259     }
260     else if (ws==2){
261     li.SetLineColor(4);
262     li.DrawLine(obt,0,obt,h1->GetMaximum());
263     }
264     }
265    
266     finalCanv->cd(1);
267     stringstream ws1, ws2;
268     ws1.str("");
269     ws2.str("");
270     ws1<<"ACQ_SETTING_MODE = 1";
271     ws2<<"ACQ_SETTING_MODE = 2";
272     TPaveText *pt=0;
273     pt = new TPaveText (.60,.92,.76,.98);
274     pt->AddText(ws1.str().c_str());
275     pt->SetTextColor(3);
276     pt->SetFillColor(10);
277     pt->SetBorderSize(0);
278     pt->Draw();
279     TPaveText *pt1=0;
280     pt1 = new TPaveText (.76,.92,.92,.98);
281     pt1->AddText(ws2.str().c_str());
282     pt1->SetTextColor(4);
283     pt1->SetFillColor(10);
284     pt1->SetBorderSize(0);
285     pt1->Draw();
286     // TPaveText *pt1=0;
287     pt1 = new TPaveText (.05,.91,.6,1);
288     pt1->AddText(oss1.str().c_str());
289     pt1->SetTextColor(1);
290     pt1->SetFillColor(10);
291     pt1->SetBorderSize(0);
292     pt1->Draw();
293 pam-rm2 1.6
294 pam-rm2 1.8 finalCanv->cd(2);
295 pam-rm2 1.6 all3->Draw();
296     all3->cd();
297     packetLength->SetStats(kFALSE);
298     packetLength->SetMarkerColor(2);
299     packetLength->SetMarkerSize(.5);
300     packetLength->SetMarkerStyle(21);
301     packetLength->GetXaxis()->SetTitle("OBT (ms)");
302     packetLength->GetXaxis()->CenterTitle();
303 pam-rm2 1.8 packetLength->GetXaxis()->SetLabelSize(0.05);
304 pam-rm2 1.6 packetLength->GetYaxis()->SetTitle("Lenght (byte)");
305     packetLength->GetYaxis()->CenterTitle();
306 pam-rm2 1.8 packetLength->GetYaxis()->SetLabelSize(0.05);
307 pam-rm2 1.6 packetLength->GetYaxis()->SetTitleSize(0.06);
308 pam-rm2 1.8 packetLength->GetYaxis()->SetTitleOffset(0.6);
309 pam-rm2 1.6 packetLength->Draw("9p");
310 pam-rm2 1.8 finalCanv->cd(2);
311     //TPaveText *pt=0;
312     pt = new TPaveText (.6,.91,.90,1);
313     pt->AddText(oss2.str().c_str());
314     pt->SetTextColor(2);
315     pt->SetFillColor(10);
316     pt->SetBorderSize(0);
317     pt->Draw();
318 pam-rm2 1.2
319 pam-rm2 1.8 finalCanv->cd(3);
320 pam-rm2 1.6 all4->Draw();
321     all4->cd();
322 pam-rm2 1.8 packeCounter->GetTitle();
323 pam-rm2 1.6 packeCounter->SetStats(kFALSE);
324     packeCounter->SetMarkerColor(4);
325 pam-rm2 1.8 packeCounter->SetMarkerSize(.2);
326 pam-rm2 1.6 packeCounter->SetMarkerStyle(21);
327     packeCounter->GetXaxis()->SetTitle("OBT (ms)");
328     packeCounter->GetXaxis()->SetTitleSize(0.05);
329     packeCounter->GetXaxis()->CenterTitle();
330 pam-rm2 1.8 packeCounter->GetXaxis()->SetLabelSize(0.05);
331 pam-rm2 1.6 packeCounter->GetYaxis()->SetTitle("Counter");
332     packeCounter->GetYaxis()->SetTitleSize(0.05);
333     packeCounter->GetYaxis()->CenterTitle();
334 pam-rm2 1.8 packeCounter->GetYaxis()->SetLabelSize(0.05);
335 pam-rm2 1.6 packeCounter->GetYaxis()->SetTitleSize(0.06);
336 pam-rm2 1.8 packeCounter->GetYaxis()->SetTitleOffset(0.6);
337     Double_t min = 500000000.;
338 pam-rm2 1.6 for (Int_t l = 0; l < nevents; l++){
339     if((packeCounter->GetBinContent(l))<= min && (packeCounter->GetBinContent(l))!= 0) min = packeCounter->GetBinContent(l);
340 pam-rm2 1.8 }
341 pam-rm2 1.6 packeCounter->SetMinimum(min);
342     packeCounter->Draw("9p");
343 pam-rm2 1.8 finalCanv->cd(3);
344     TPaveText *pt2=0;
345     pt2 = new TPaveText (.6,.91,.90,1);
346     pt2->AddText(oss3.str().c_str());
347     pt2->SetTextColor(4);
348     pt2->SetFillColor(10);
349     pt2->SetBorderSize(0);
350     pt2->Draw();
351 pam-rm2 1.4
352    
353 pam-rm2 1.6 //-------------------------------------------------
354     TList *list = new TList;
355     Int_t numkey;
356     TObject *key = new TObject;
357     const char *name;
358     TTree* tr = new TTree;
359     Long64_t nevntskey=0;
360     list = file->GetListOfKeys();
361     numkey = file->GetNkeys();
362 pam-rm2 1.8 int ll=0;
363     int kk=0;
364    
365 pam-rm2 1.6 for (Int_t i=0; i<numkey; i++){
366     key = list->At(i);
367     name=(char *)(key->GetName());
368     tr = (TTree*)file->Get(name);
369     if (tr->IsZombie()) continue;
370     tr->SetBranchAddress("Header", &eh);
371     TBranch *Br = tr->GetBranch("Header");
372     nevntskey = tr->GetEntries();
373    
374 pam-rm2 1.8 if(nevntskey !=0){
375     Int_t size1=nevntskey;
376     Double_t *PscuCounter1 = new Double_t[size1];
377     Double_t *FileOffset1 = new Double_t[size1];
378     Double_t *tempo1 = new Double_t[size1];
379     //Double_t *salto= new Double_t[100];
380     //Double_t *salto1= new Double_t[100];
381    
382     int l=0;
383     while(l<nevntskey){
384     Br->GetEntry(l);
385     ph = eh->GetPscuHeader();
386     PscuCounter1[l]= ph->GetCounter();
387     FileOffset1[l]=ph->GetFileOffset();
388     tempo1[l]=ph->GetOrbitalTime();
389     if(ph->GetFileOffset()<= primoffset){
390     primoffset =ph->GetFileOffset();
391     primotempo=ph->GetOrbitalTime();
392     }
393     if(ph->GetFileOffset()>= ultimoffset){
394     ultimoffset = ph->GetFileOffset();
395     ultimotempo = ph->GetOrbitalTime();
396     }
397     /*if(l>0){
398     if(tempo1[l] < tempo1[l-1]){
399     salto1[kk] =(const Double_t*)tempo1[l-1];
400     salto[ll]= ph->GetOrbitalTime();
401     ll=ll+1;
402     kk=kk+1;
403 pam-rm2 1.6 }
404 pam-rm2 1.8 }*/
405     l++;
406     }
407    
408     TGraph *graph3= new TGraph(nevntskey, (const Double_t*)FileOffset1, (const Double_t*)PscuCounter1);
409     graph3->SetMarkerColor(3);
410     graph3->SetMarkerSize(.2);
411     graph3->SetMarkerStyle(21);
412     mg1->Add(graph3);
413    
414     TGraph *graph4= new TGraph(nevntskey, (const Double_t*)FileOffset1, (const Double_t*)tempo1);
415     graph4->SetMarkerColor(kBlue);
416     graph4->SetMarkerSize(.2);
417     graph4->SetMarkerStyle(21);
418     mg2->Add(graph4);
419 pam-rm2 1.6 }
420     }
421 pam-rm2 1.8
422     kk=kk-1;
423     ll=ll-1;
424    
425     TLatex *lat = new TLatex();
426     lat->SetTextAlign(12);
427     lat->SetTextSize(0.15);
428     lat->SetTextColor(kBlue);
429     //------------ Second Canvas ---------------------//
430     finalCanv->cd(4);
431 pam-rm2 1.6 all1->Draw();
432     all1->cd();
433 pam-rm2 1.8
434     oss1.str("");
435     oss1 << filename.Data() <<": PscuCounter vs FileOffset.";
436     //mg1->SetTitle(oss1.str().c_str());
437 pam-rm2 1.3 mg1->Draw("AP");
438     mg1->GetXaxis()->SetTitle("File Offset");
439     mg1->GetXaxis()->CenterTitle();
440 pam-rm2 1.8 mg1->GetXaxis()->SetTitleSize(0.045);
441     mg1->GetXaxis()->SetLabelSize(0.05);
442 pam-rm2 1.3 mg1->GetYaxis()->SetTitle("Counter");
443     mg1->GetYaxis()->CenterTitle();
444 pam-rm2 1.4 mg1->GetYaxis()->SetTitleSize(0.05);
445 pam-rm2 1.8 mg1->GetYaxis()->SetLabelSize(0.06);
446     mg1->GetYaxis()->SetTitleOffset(0.7);
447     finalCanv->cd(4);
448     TPaveText *pt3=0;
449     pt3 = new TPaveText (.60,.91,.90,1);
450     pt3->AddText(oss1.str().c_str());
451     pt3->SetTextColor(3);
452     pt3->SetFillColor(10);
453     pt3->SetBorderSize(0);
454     pt3->Draw();
455    
456 pam-rm2 1.4
457 pam-rm2 1.8 finalCanv->cd(5);
458 pam-rm2 1.6 all->Draw();
459     all->cd();
460    
461     oss3.str("");
462 pam-rm2 1.8 oss3 << filename.Data() <<" OBT vs FileOffset. First packet = "<<primotempo <<" ms, Last packet = "<<ultimotempo<<" ms.";
463     //mg2->SetTitle(oss3.str().c_str());
464 pam-rm2 1.4 mg2->Draw("AP");
465     mg2->GetXaxis()->SetTitle("File Offset");
466     mg2->GetXaxis()->CenterTitle();
467 pam-rm2 1.8 mg2->GetXaxis()->SetTitleSize(0.045);
468     mg2->GetXaxis()->SetLabelSize(0.05);
469 pam-rm2 1.4 mg2->GetYaxis()->SetTitle("OBT");
470     mg2->GetYaxis()->CenterTitle();
471     mg2->GetYaxis()->SetTitleSize(0.05);
472 pam-rm2 1.8 mg2->GetYaxis()->SetLabelSize(0.055);
473     mg2->GetYaxis()->SetTitleOffset(0.6);
474     /*pad->Draw();
475     pad->cd();
476     stringa.str("");
477     stringa << "new data from "<< primotempo;//<<" to "<< salto1;
478     lat->DrawLatex(0.08, 0.8, stringa.str().c_str());
479 pam-rm2 1.6 double min2 = 0.8;
480 pam-rm2 1.8 if(ll > 0){
481 pam-rm2 1.6 pad->Draw();
482     pad->cd();
483 pam-rm2 1.8 for(Int_t k=0; k <(kk); k++){
484 pam-rm2 1.6 stringa.str("");
485 pam-rm2 1.8 stringa << "old data from "<< salto[kk]<< " to "<<salto1[kk+1];
486 pam-rm2 1.6 min2=min2-0.1;
487     lat->DrawLatex(0.08, min2,stringa.str().c_str());
488     //cout<<salto[kk]<<";\n";
489     }
490 pam-rm2 1.8 }*/
491 pam-rm2 1.6
492 pam-rm2 1.8 finalCanv->cd(5);
493     TPaveText *pt4=0;
494     pt4 = new TPaveText (.38 ,.91,.92,1);
495     pt4->AddText(oss3.str().c_str());
496     pt4->SetTextColor(kBlue);
497     pt4->SetFillColor(10);
498     pt4->SetBorderSize(0);
499     pt4->Draw();
500 pam-rm2 1.6
501 pam-rm2 1.8 // finalCanv->Draw();
502     finalCanv->Update();
503 pam-rm2 1.6
504 pam-rm2 1.2 oss1.str("");
505     oss1 << outDir.Data() << filename.Data();
506 pam-rm2 1.8 oss1 << "_HeaderScan"<<"." << format.Data();
507    
508     finalCanv->SaveAs(oss1.str().c_str());
509 pam-rm2 1.2
510 pam-rm2 1.4
511 pam-rm2 1.8
512 pam-rm2 1.3 file->Close();
513 pam-rm2 1.1
514     }
515    
516     int main(int argc, char* argv[]){
517     TString path;
518     TString outDir = "./";
519     TString format = "jpg";
520     if (argc < 2){
521     printf("You have to insert at least the file to analyze \n");
522     printf("Try '--help' for more information. \n");
523     exit(1);
524     }
525     if (!strcmp(argv[1], "--help")){
526     printf( "Usage: HeaderScan FILE [OPTION] \n");
527     printf( "\t --help Print this help and exit \n");
528     printf( "\t -outDir[path] Path where to put the output [default ./] \n");
529     printf( "\t -format[jpg|ps|gif] Format for output files [default 'jpg'] \n");
530     exit(1);
531     }
532     path=argv[1];
533 pam-rm2 1.2 for (int i = 2; i < argc; i++){
534 pam-rm2 1.1 if (!strcmp(argv[i], "-outDir")){
535     if (++i >= argc){
536     printf( "-outDir needs arguments. \n");
537     printf( "Try '--help' for more information. \n");
538     exit(1);
539     }
540     else{
541     outDir = argv[i];
542     continue;
543     }
544 pam-rm2 1.2 }
545 pam-rm2 1.1 if (!strcmp(argv[i], "-format")){
546     if (++i >= argc){
547     printf( "-format needs arguments. \n");
548     printf( "Try '--help' for more information. \n");
549     exit(1);
550     }
551     else{
552     format = argv[i];
553     continue;
554     }
555     }
556     }
557     HeaderScan(argv[1], outDir, format);
558     }
559 pam-rm2 1.7

  ViewVC Help
Powered by ViewVC 1.1.23