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

Contents of /quicklook/QLflightTmtc_Header/HeaderScan.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.12 - (show annotations) (download)
Tue Oct 24 08:23:17 2006 UTC (18 years, 1 month ago) by pam-rm2
Branch: MAIN
CVS Tags: v3r03
Changes since 1.11: +0 -1 lines
aggiornato il PacketScan perche' sia compatibile con Yoda_3/14

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

  ViewVC Help
Powered by ViewVC 1.1.23