/[PAMELA software]/tof/flight/ToFReprocessing/ToFreproc.cpp
ViewVC logotype

Annotation of /tof/flight/ToFReprocessing/ToFreproc.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.2 - (hide annotations) (download)
Thu Feb 5 09:12:35 2009 UTC (15 years, 9 months ago) by mocchiut
Branch: MAIN
Changes since 1.1: +191 -186 lines
Small changes

1 mocchiut 1.1 #if !defined(__CINT__) || defined(__MAKECINT__)
2    
3    
4     #include <fstream>
5     #include <sstream>
6     #include <string>
7     #include <math.h>
8     #include <iostream>
9     #include <TTree.h>
10     #include <TClassEdit.h>
11     #include <TObject.h>
12     #include <TList.h>
13     #include <TSystem.h>
14     #include <TSystemDirectory.h>
15     #include <TString.h>
16     #include <TFile.h>
17     #include <TClass.h>
18     #include <TCanvas.h>
19     #include <TH1.h>
20     #include <TH1F.h>
21     #include <TH2D.h>
22     #include <TLatex.h>
23     #include <TPad.h>
24     #include <TPaveLabel.h>
25     #include <TChain.h>
26     #include <TGraph.h>
27     #include <TMultiGraph.h>
28     #include <TLine.h>
29     #include <TStyle.h>
30     #include <TF1.h>
31     #include <TH1F.h>
32     #include <TTree.h>
33     #include <TCanvas.h>
34     #include <PamLevel2.h>
35     #include <iomanip>
36     //
37     #include <ToFCore.h>
38     //
39     // Declaration of the core fortran routines
40     //
41     #define tofl2com tofl2com_
42     extern "C" int tofl2com();
43     #define toftrk toftrk_
44     extern "C" int toftrk();
45     ///extern "C" int toftrk(float);
46     #define rdtofcal rdtofcal_
47     extern "C" int rdtofcal(char [], int *);
48    
49     using namespace std;
50    
51     #endif
52    
53     //===================================
54     // global variables
55     //===================================
56     TString DIR;
57     TString OUTF;
58     TString LIST;
59     TString OPTIONS;
60    
61     PamLevel2 *pam_event = NULL;
62    
63     //==========================================
64     //000000000000000000000000000000000000000000
65     //==========================================
66     Int_t Loop(TString ddir,TString list, TString options){
67     //
68     //
69     extern struct ToFInput tofinput_;
70     extern struct ToFOutput tofoutput_;
71     // --------------------
72     // read input file/list
73     // --------------------
74     pam_event = new PamLevel2(ddir,list,options);
75 mocchiut 1.2 pam_event->SetSELLI(2);
76 mocchiut 1.1 TTree::SetMaxTreeSize(1000*Long64_t(2000000000));
77 mocchiut 1.2 //
78 mocchiut 1.1 TString outfile_t =0;
79     outfile_t = OUTF;
80     outfile_t.Append("b.root");
81     TFile *outft = (TFile*)gROOT->FindObject(outfile_t);
82     if (outft) outft->Close();
83     outft = new TFile(outfile_t,"RECREATE");
84     if(outft->IsZombie()){
85     cout << "Output file could not be created\n";
86     return 1;
87     };
88     cout << "Created output file: "<<outft->GetName()<<endl;
89     //
90     outft->cd();
91     ToFLevel2 *tof = new ToFLevel2();
92     TTree *toft = new TTree("ToF","PAMELA Level2 ToF data");
93     toft->SetAutoSave(900000000000000LL);
94     tof->Set();//ELENA **TEMPORANEO?**
95     toft->Branch("ToFLevel2","ToFLevel2",&tof);
96     //
97     pam_event->CreateCloneTrees(outft);
98     //
99     Int_t ntrkentry = 0;
100     Int_t npmtentry = 0;
101     //
102     ULong64_t nevents = pam_event->GetEntries();
103     printf("\n\n Running on %llu events \n\n",nevents);
104     //
105     TFile *tfile = TFile::Open(list.Data(),"READ");
106     TTree *toftr = (TTree*)tfile->Get("ToF");
107     ToFLevel2 *tofl2 = new ToFLevel2();
108     toftr->SetBranchAddress("ToFLevel2", &tofl2);
109     ULong64_t ntofev = toftr->GetEntries();
110     printf(" ToF events are %llu \n",ntofev);
111     //
112     //
113     GL_PARAM *glparam = new GL_PARAM();
114     TrkLevel2 *trk = new TrkLevel2();
115     //
116     TString host;
117     TString user;
118     TString psw;
119     const char *pamdbhost=gSystem->Getenv("PAM_DBHOST");
120     const char *pamdbuser=gSystem->Getenv("PAM_DBUSER");
121     const char *pamdbpsw=gSystem->Getenv("PAM_DBPSW");
122     if ( !pamdbhost ) pamdbhost = "";
123     if ( !pamdbuser ) pamdbuser = "";
124     if ( !pamdbpsw ) pamdbpsw = "";
125     if ( strcmp(pamdbhost,"") ) host = pamdbhost;
126     if ( strcmp(pamdbuser,"") ) user = pamdbuser;
127     if ( strcmp(pamdbpsw,"") ) psw = pamdbpsw;
128     //
129     //
130     TSQLServer *dbc = TSQLServer::Connect(host.Data(),user.Data(),psw.Data());
131     if ( !dbc->IsConnected() ) return 1;
132     stringstream myquery;
133     myquery.str("");
134     myquery << "SET time_zone='+0:00'";
135     dbc->Query(myquery.str().c_str());
136     glparam->Query_GL_PARAM(1,1,dbc); // parameters stored in DB in GL_PRAM table
137     trk->LoadField(glparam->PATH+glparam->NAME);
138     //
139     Long64_t myrun = 0;
140     Long64_t mysrun = -1;
141     Int_t adc[4][12];
142     Int_t tdc[4][12];
143     Float_t tdcc[4][12];
144     //
145     Bool_t defcal = true;
146     for(ULong64_t iev=0; iev<nevents; iev++){
147    
148    
149 mocchiut 1.2 //==================================================================
150 mocchiut 1.1
151    
152     // for(ULong64_t iev=0; iev<50; iev++){
153     //
154     if ( !(iev%1000) ) printf(" %iK \n",(int)iev/1000);
155     pam_event->Clear();
156     //
157     if( pam_event->GetEntry(iev) ){
158     //
159     pam_event->FillCloneTrees();
160     toftr->GetEntry(iev);
161     //
162     // variable to save information about the tof calibration used
163     //
164     //
165     myrun=pam_event->GetRunInfo()->ID;
166     if ( myrun != mysrun ){
167     // printf(" rhtime %u myrun %i mysrun %i \n",pam_event->GetRunInfo()->RUNHEADER_TIME,(int)myrun,(int)mysrun);
168     mysrun = myrun;
169     //
170     Int_t error=glparam->Query_GL_PARAM(pam_event->GetRunInfo()->RUNHEADER_TIME,201,dbc); // parameters stored in DB in GL_PRAM table
171     if ( error<0 ) {
172     return(1);
173     };
174     //
175     printf(" Reading ToF parameter file: %s \n",(glparam->PATH+glparam->NAME).Data());
176     //
177     if ( (UInt_t)glparam->TO_TIME != (UInt_t)4294967295UL ) defcal = false;
178     //
179     Int_t nlen = (Int_t)(glparam->PATH+glparam->NAME).Length();
180     rdtofcal((char *)(glparam->PATH+glparam->NAME).Data(),&nlen);
181     //
182     };
183     //
184    
185    
186 mocchiut 1.2 //================================================================
187     //==================================================================
188     //--- Define absolute time
189     Float_t tabs=pam_event->GetOrbitalInfo()->absTime;
190    
191     //==================================================================
192     Float_t dedx_corr_m[2000][48],dedx_corr[48];
193     Double_t mtime[2000],t1,t2,tm;
194     Float_t yhelp1,yhelp2,slope,inter,thelp1,thelp2;
195     Float_t xmean1,xwidth1;
196    
197     Int_t ical,ii,j,jj;
198    
199     if (iev==0) {
200    
201     ical=0; // counter set to zero if first-time reading
202    
203     //-----------------------------------------------------------
204     // Here I read the dEdx_korr parameters
205     //-----------------------------------------------------------
206    
207     jj=0;
208    
209     ifstream fin("adcmonitor.7th.100k.cut.dat");
210    
211     while (! fin.eof()) {
212     fin>>t1>>tm>>t2;
213     // cout<<jj<<" "<<tm<<endl;
214     cout << setiosflags(ios::fixed) << setw(10) << setprecision(3) << tm << endl;
215     mtime[jj]=tm;
216     for (ii=0; ii<48;ii++) {
217     fin>>j>>xmean1>>xwidth1;
218     dedx_corr_m[jj][ii]=xmean1;
219     }
220     jj=jj+1;
221     // printf(" kk %i \n",jj);
222     }
223     // printf(" 1ical %i \n",ical);
224    
225     fin.close();
226    
227    
228     while (tabs<mtime[ical] || tabs > mtime[ical+1]) {
229     // printf(" ical %i \n",ical);
230     ical = ical+1;
231     }
232     // ical = ical-1;
233     cout<<"abs time "<<tabs<<" limit low "<<mtime[ical]<<" limit up "<<mtime[ical+1]<<" ical "<<ical<<endl;
234     } // if iev==0...
235    
236     //==================================================================
237     //== End of first time reading & filling the array
238     //==================================================================
239    
240     //==================================================================
241     //== if time is outside time limits:
242     //==================================================================
243    
244     if (tabs<mtime[ical] || tabs>mtime[ical+1]) {
245     cout<<"Checking Time Limits!"<<endl;
246     ical=0;
247     while (tabs > mtime[ical+1] || tabs<mtime[ical]) {
248     ical = ical+1;
249     }
250     // ical = ical-1;
251     cout<<"abs time "<<tabs<<" limit low "<<mtime[ical]<<" limit up "<<mtime[ical+1]<<" ical "<<ical<<endl;
252     };
253 mocchiut 1.1
254 mocchiut 1.2 // printf(" 2ical %i \n",ical);
255     //==================================================================
256     //== interpolate betwen time limits
257     //==================================================================
258    
259     thelp1 = mtime[ical];
260     thelp2 = mtime[ical+1];
261    
262     for (ii=0; ii<48;ii++) {
263     yhelp1 = dedx_corr_m[ical][ii];
264     yhelp2 = dedx_corr_m[ical+1][ii];
265     slope = (yhelp2-yhelp1)/(thelp2-thelp1);
266     inter = yhelp1 - slope*thelp1;
267     dedx_corr[ii] = slope*tabs + inter;
268     // if (ii==0) cout<<thelp1<<" "<<thelp2<<" "<<tabs<<" "<<yhelp1<<" "<<yhelp2<<" "<<dedx_corr[0]<<endl;
269     }
270 mocchiut 1.1
271 mocchiut 1.2 //================================================================
272     //================================================================
273     // printf("cpippo \n");
274 mocchiut 1.1
275     // process tof data
276     //
277     for (Int_t hh=0; hh<12;hh++){
278     for (Int_t kk=0; kk<4;kk++){
279 mocchiut 1.2 adc[kk][hh] = 4095;
280     tdc[kk][hh] = 4095;
281     tdcc[kk][hh] = 4095.;
282     tofinput_.adc[hh][kk] = 4095;
283     tofinput_.tdc[hh][kk] = 4095;
284 mocchiut 1.1 };
285     };
286 mocchiut 1.2 // memset(adc, 0, 12*4*sizeof(Int_t));
287     // memset(tdc, 0, 12*4*sizeof(Int_t));
288 mocchiut 1.1 Int_t gg = 0;
289     Int_t hh = 0;
290     Int_t adcf[48];
291     memset(adcf, 0, 48*sizeof(Int_t));
292     Int_t tdcf[48];
293     memset(tdcf, 0, 48*sizeof(Int_t));
294     for (Int_t pm=0; pm < tofl2->ntrk() ; pm++){
295     ToFTrkVar *ttf = tofl2->GetToFTrkVar(pm);
296     for ( Int_t nc=0; nc < ttf->npmttdc; nc++){
297 mocchiut 1.2 if ( (ttf->tdcflag).At(nc) != 0 ) tdcf[(ttf->pmttdc).At(nc)] = 1;
298 mocchiut 1.1 };
299     for ( Int_t nc=0; nc < ttf->npmtadc; nc++){
300 mocchiut 1.2 if ( (ttf->adcflag).At(nc) != 0 ) adcf[(ttf->pmtadc).At(nc)] = 1;
301 mocchiut 1.1 };
302     };
303     //
304     for (Int_t pm=0; pm < tofl2->npmt() ; pm++){
305     ToFPMT *pmt = tofl2->GetToFPMT(pm);
306     tofl2->GetPMTIndex(pmt->pmt_id, gg, hh);
307     if ( adcf[pmt->pmt_id] == 0 ){
308 mocchiut 1.2 tofinput_.adc[gg][hh] = (int)pmt->adc;
309     adc[hh][gg] = (int)pmt->adc;
310 mocchiut 1.1 };
311     if ( tdcf[pmt->pmt_id] == 0 ){
312 mocchiut 1.2 tofinput_.tdc[gg][hh] = (int)pmt->tdc;
313     tdc[hh][gg] = (int)pmt->tdc;
314 mocchiut 1.1 };
315     tdcc[hh][gg] = (float)pmt->tdc_tw;
316     // Int_t pppid = tof->GetPMTid(hh,gg);
317     // printf(" pm %i pmt_id %i pppid %i hh %i gg %i tdcc %f tdc %f adc %f \n",pm,pmt->pmt_id,pppid,hh,gg,pmt->tdc_tw,pmt->tdc,pmt->adc);
318     };
319    
320    
321 mocchiut 1.2 // for (Int_t pm=0; pm <48 ; pm++){
322     // tof->GetPMTIndex(pm, gg, hh);
323     // tofinput_.tdc[hh][gg] = (int)500.;
324     // tofinput_.adc[hh][gg] = (int)500.;
325     // tdc[hh][gg] = (int)500.;
326     // adc[hh][gg] = (int)500.;
327     // // printf(" hh %i gg %i tdc %f adc %f \n",hh,gg,pmt->tdc,pmt->adc);
328     // };
329 mocchiut 1.1 //
330     for (Int_t hh=0; hh<5;hh++){
331     tofinput_.patterntrig[hh]=pam_event->GetTrigLevel2()->patterntrig[hh];
332     };
333     //
334     tof->Clear();
335     Int_t pmt_id = 0;
336     ToFPMT *t_pmt = new ToFPMT();
337     if(!(tof->PMT))tof->PMT = new TClonesArray("ToFPMT",12); //ELENA
338     TClonesArray &tpmt = *tof->PMT;
339     ToFTrkVar *t_tof = new ToFTrkVar();
340     if(!(tof->ToFTrk))tof->ToFTrk = new TClonesArray("ToFTrkVar",2); //ELENA
341     TClonesArray &t = *tof->ToFTrk;
342     //
343     //
344     // Here we have calibrated data, ready to be passed to the FORTRAN routine which will extract common and track-related variables.
345     //
346     npmtentry = 0;
347     //
348     ntrkentry = 0;
349     //
350     // Calculate tracks informations from ToF alone
351     //
352     tofl2com();
353     //
354     memcpy(tof->tof_j_flag,tofoutput_.tof_j_flag,6*sizeof(Int_t));
355     //
356     t_tof->trkseqno = -1;
357     //
358     // and now we must copy from the output structure to the level2 class:
359     //
360     t_tof->npmttdc = 0;
361     //
362     for (Int_t hh=0; hh<12;hh++){
363     for (Int_t kk=0; kk<4;kk++){
364     if ( tofoutput_.tofmask[hh][kk] != 0 ){
365     pmt_id = tof->GetPMTid(kk,hh);
366     t_tof->pmttdc.AddAt(pmt_id,t_tof->npmttdc);
367     t_tof->tdcflag.AddAt(tofoutput_.tdcflagtof[hh][kk],t_tof->npmttdc); // gf: Jan 09/07
368     t_tof->npmttdc++;
369     };
370     };
371     };
372     for (Int_t kk=0; kk<13;kk++){
373     t_tof->beta[kk] = tofoutput_.betatof_a[kk];
374     }
375     //
376     t_tof->npmtadc = 0;
377     for (Int_t hh=0; hh<12;hh++){
378     for (Int_t kk=0; kk<4;kk++){
379     if ( tofoutput_.adctof_c[hh][kk] < 1000 ){
380     // t_tof->dedx.AddAt(tofoutput_.adctof_c[hh][kk],t_tof->npmtadc); // EMILIANO
381     pmt_id = tof->GetPMTid(kk,hh);
382     t_tof->dedx.AddAt((tofoutput_.adctof_c[hh][kk]*4./dedx_corr[pmt_id]),t_tof->npmtadc); // EMILIANO
383     t_tof->pmtadc.AddAt(pmt_id,t_tof->npmtadc);
384     t_tof->adcflag.AddAt(tofoutput_.adcflagtof[hh][kk],t_tof->npmtadc); // gf: Jan 09/07
385     t_tof->npmtadc++;
386     };
387     };
388     };
389     //
390     memcpy(t_tof->xtofpos,tofoutput_.xtofpos,sizeof(t_tof->xtofpos));
391     memcpy(t_tof->ytofpos,tofoutput_.ytofpos,sizeof(t_tof->ytofpos));
392     memcpy(t_tof->xtr_tof,tofoutput_.xtr_tof,sizeof(t_tof->xtr_tof));
393     memcpy(t_tof->ytr_tof,tofoutput_.ytr_tof,sizeof(t_tof->ytr_tof));
394     //
395     new(t[ntrkentry]) ToFTrkVar(*t_tof);
396     ntrkentry++;
397     t_tof->Clear();
398     //
399     //
400     //
401     t_pmt->Clear();
402     //
403     for (Int_t hh=0; hh<12;hh++){
404     for (Int_t kk=0; kk<4;kk++){
405 mocchiut 1.2 // new WM
406     // if ( tofoutput_.tdc_c[hh][kk] < 4095 || adc[kk][hh] < 4095 || tdc[kk][hh] < 4095 ){
407 mocchiut 1.1 if ( tdcc[kk][hh] < 4095. || adc[kk][hh] < 4095 || tdc[kk][hh] < 4095 ){
408     //
409     t_pmt->pmt_id = tof->GetPMTid(kk,hh);
410     t_pmt->tdc_tw = tofoutput_.tdc_c[hh][kk];
411     t_pmt->adc = (Float_t)adc[kk][hh];
412     t_pmt->tdc = (Float_t)tdc[kk][hh];
413     //
414     new(tpmt[npmtentry]) ToFPMT(*t_pmt);
415     npmtentry++;
416     t_pmt->Clear();
417     };
418     };
419     };
420     //
421     // Calculate track-related variables
422     //
423     if ( pam_event->GetTrkLevel2()->ntrk() > 0 ){
424     //
425     // We have at least one track
426     //
427     //
428     // Run over tracks
429     //
430     for(Int_t nt=0; nt < pam_event->GetTrkLevel2()->ntrk(); nt++){
431     //
432     TrkTrack *ptt = pam_event->GetTrkLevel2()->GetStoredTrack(nt);
433     //
434     // Copy the alpha vector in the input structure
435     //
436     for (Int_t e = 0; e < 5 ; e++){
437     tofinput_.al_pp[e] = ptt->al[e];
438     };
439     //
440     // Get tracker related variables for this track
441     //
442     toftrk();
443 mocchiut 1.2 // toftrk(thelp);
444 mocchiut 1.1 //
445     // Copy values in the class from the structure (we need to use a temporary class to store variables).
446     //
447     t_tof->npmttdc = 0;
448     for (Int_t hh=0; hh<12;hh++){
449     for (Int_t kk=0; kk<4;kk++){
450     if ( tofoutput_.tofmask[hh][kk] != 0 ){
451     pmt_id = tof->GetPMTid(kk,hh);
452     t_tof->pmttdc.AddAt(pmt_id,t_tof->npmttdc);
453     t_tof->tdcflag.AddAt(tofoutput_.tdcflag[hh][kk],t_tof->npmttdc); // gf: Jan 09/07
454     t_tof->npmttdc++;
455     };
456     };
457     };
458     for (Int_t kk=0; kk<13;kk++){
459     t_tof->beta[kk] = tofoutput_.beta_a[kk];
460     };
461     //
462     t_tof->npmtadc = 0;
463     for (Int_t hh=0; hh<12;hh++){
464     for (Int_t kk=0; kk<4;kk++){
465     if ( tofoutput_.adc_c[hh][kk] < 1000 ){
466 mocchiut 1.2 // t_tof->dedx.AddAt(tofoutput_.adc_c[hh][kk],t_tof->npmtadc); // EMILIANO
467 mocchiut 1.1 pmt_id = tof->GetPMTid(kk,hh);
468     t_tof->dedx.AddAt((tofoutput_.adc_c[hh][kk]*4./dedx_corr[pmt_id]),t_tof->npmtadc); // EMILIANO
469     t_tof->pmtadc.AddAt(pmt_id,t_tof->npmtadc);
470     t_tof->adcflag.AddAt(tofoutput_.adcflag[hh][kk],t_tof->npmtadc); // gf: Jan 09/07
471     t_tof->npmtadc++;
472     };
473     };
474     };
475     //
476     memcpy(t_tof->xtofpos,tofoutput_.xtofpos,sizeof(t_tof->xtofpos));
477     memcpy(t_tof->ytofpos,tofoutput_.ytofpos,sizeof(t_tof->ytofpos));
478     memcpy(t_tof->xtr_tof,tofoutput_.xtr_tof,sizeof(t_tof->xtr_tof));
479     memcpy(t_tof->ytr_tof,tofoutput_.ytr_tof,sizeof(t_tof->ytr_tof));
480     //
481     // Store the tracker track number in order to be sure to have shyncronized data during analysis
482     //
483     t_tof->trkseqno = nt;
484     //
485     // create a new object for this event with track-related variables
486     //
487     new(t[ntrkentry]) ToFTrkVar(*t_tof);
488     ntrkentry++;
489     t_tof->Clear();
490     //
491     }; // loop on all the tracks
492     };
493     //
494     tof->unpackError = tofl2->unpackError;
495     if ( defcal ){
496     tof->default_calib = 1;
497     } else {
498     tof->default_calib = 0;
499     };
500     //
501     // Fill the rootple
502     //
503     toft->Fill();
504     //
505     // toft->Show();
506     // toftr->Show();
507     //
508     delete t_tof;
509     //
510     //
511     };
512     //
513     };
514     //
515     pam_event->Clear();
516     //
517     outft->cd();
518     pam_event->WriteCloneTrees();
519     toftr->Delete();
520    
521     toft->Write();
522     outft->Close();
523     //
524     return 0;
525     //
526     }
527    
528     /////////////////////////////////////////////////////////////////
529    
530     #if !defined(__CINT__)
531    
532     void usage(){
533    
534 mocchiut 1.2 cout << "------------------------------------------------------------"<<endl;
535     cout << "Loop over events (on one or more Level2-files), applying some selection cuts (defined in My-Selection.cpp), \n";
536     cout << "creating output histograms (defined in My-Histos.cpp) and/or trees with selected events. \n \n ";
537     cout << "USAGE:"<<endl;
538     cout << "-processDir DIR - Level2 data directory \n";
539     cout << "-processList LIST - list of files (.txt) or single file (.root) to be analysed \n";
540     cout << "-outputFile PATH - name of the output file \n";
541     cout << "-NumEvents XXX - number of events to be analysed \n";
542     cout << "--debug, -g - debug mode \n";
543     cout << "--help, -h - print this help \n";
544     cout << "-options [ options ] - options: \n";
545     cout << " fillHistos --> create an output file with histograms \n";
546     cout << " fillTree --> create an output file with trees storing the selected events \n ";
547     cout << " +(-)ALL --> inlcude(exclude) all trees and branches \n " ;
548     cout << " +(-)TRK1 +(-)TRK2 +(-)CAL1 +(-)CAL2 +(-)TOF +(-)TRG +(-)ND +(-)S4 +(-)ORB --> inlcude(exclude) trees and branches \n" ;
549     cout << "------------------------------------------------------------"<<endl;
550 mocchiut 1.1 }
551     //
552     int HandleInputPar(int argc, char **argv){
553    
554 mocchiut 1.2 if(argc>1){
555 mocchiut 1.1
556 mocchiut 1.2 if(!strcmp(argv[1], "-h") || !strcmp(argv[1], "--help") ){
557     usage();
558     return(1);
559     };
560     // -----------------------
561     // Read input parameters
562     // -----------------------
563     DIR = gSystem->WorkingDirectory();
564     LIST = "";
565     OUTF = "myfile";
566     OPTIONS = "+AUTO -TOF";
567 mocchiut 1.1
568 mocchiut 1.2 for (int i = 1; i < argc; i++){
569     // -----------------------------------------------------//
570     if (!strcmp(argv[i], "-processDir")){
571     if (++i >= argc) throw -1;
572     DIR = argv[i];
573     cout << "processDir "<<DIR<<endl;
574     continue;
575     }
576     // -----------------------------------------------------//
577     else if (!strcmp(argv[i], "-processList")){
578     if (++i >= argc) throw -1;
579     LIST = argv[i];
580     cout << "processList "<<LIST<<endl;
581     continue;
582     }
583     // -----------------------------------------------------//
584     else if (!strcmp(argv[i], "-outputFile")){
585     if (++i >= argc) throw -1;
586     OUTF = argv[i];
587     cout << "outputFile "<<OUTF<<endl;
588     continue;
589     }
590     // -----------------------------------------------------//
591     else if (!strcmp(argv[i], "-options")){
592     if (++i >= argc) throw -1;
593     OPTIONS = argv[i];
594     if( OPTIONS.Contains("[") ){
595     do{
596     if (++i >= argc) throw -1;
597     OPTIONS.Append(argv[i]);
598     }while(!OPTIONS.Contains("]"));
599     }else cout << "wrong option format --> ignoring " << endl;
600     }
601     else{
602     cout << "Unidentified input parameter. Ignored."<< endl;
603     };
604 mocchiut 1.1 };
605 mocchiut 1.2 }else{
606     usage();
607     return(1);
608     };
609     // -----------------------
610     // Check input parameters
611     // -----------------------
612 mocchiut 1.1
613    
614 mocchiut 1.2 return(0);
615 mocchiut 1.1
616     };
617     //
618    
619     int main(int argc, char **argv)
620     {
621    
622 mocchiut 1.2 if( HandleInputPar(argc,argv) )return(1);
623 mocchiut 1.1
624 mocchiut 1.2 cout << "OPTIONS "<<OPTIONS<<endl;
625     Loop(DIR,LIST,OPTIONS);
626 mocchiut 1.1
627 mocchiut 1.2 cout << "Back to main - end"<<endl;
628 mocchiut 1.1
629 mocchiut 1.2 return 0;
630 mocchiut 1.1
631     };
632    
633     #endif

  ViewVC Help
Powered by ViewVC 1.1.23