/[PAMELA software]/quicklook/QLflightS4_ND/S4_QL.cpp
ViewVC logotype

Contents of /quicklook/QLflightS4_ND/S4_QL.cpp

Parent Directory Parent Directory | Revision Log Revision Log


Revision 1.2 - (show annotations) (download)
Fri Jun 23 16:06:09 2006 UTC (18 years, 5 months ago) by pam-rm2
Branch: MAIN
CVS Tags: v1r02, v1r03
Changes since 1.1: +24 -15 lines
corretto il bakground del ND e la rate media di S4

1
2 /*
3 * S4 Quick Look
4 * Author Marcelli-Malvezzi
5 * Version 1.00 - March 2006
6 *
7 * Description - The aim of S4 QL software is to monitor the bahaviour of this detector.
8 * It creates three canvases: the first one contains the histogram of S4 data, the second one
9 * is relative to the time behaviour of collected data while the last shows S4 rate
10 * (information from trigger Packet).
11 * See documentation for a more detailed description of the output.
12 *
13 *
14 * Parameters:
15 * TSTring base - the path to the root directory for the specific Pamela unpack session
16 * There is no default value, without this input the program will not run
17 * TString outDir - the path where to save the output image (Default = ./)
18 * TString format - the format which will be used fo rsave the produced images (Default = "jpg")
19 * Float_t DeltaT - the time interval in minute for calculation of average S4 data for minute,
20 * see S4_QL_2 plot (Default = 1 minute)
21 *
22 *
23 * Version 1.1 - June 2006
24 * Fixed bugs: the vector "trcss" was inizilized to have dimension 10, but for files with large number of events this
25 * is not true; now this vector is inizialized at 100.
26 *
27 * the threshold at which S4 is set and the trigger configuration can change in the file; all these changes are reported
28 * in a pad
29 *
30 * for a large namber of events is not possible to have vectors, so all graphs have been converted in histograms
31 *******/
32
33
34 #include <iostream>
35 #include <fstream>
36 #include <sstream>
37 #include <math.h>
38 #include "TStyle.h"
39 #include "TFile.h"
40 #include "TList.h"
41 #include "TTree.h"
42 #include "TLatex.h"
43 #include "TObjString.h"
44 #include "TCanvas.h"
45 #include "TGraph.h"
46 #include "TH1F.h"
47 #include "TF1.h"
48 #include "TGaxis.h"
49 #include "TString.h"
50 #include "TPaveText.h"
51 #include "EventHeader.h"
52 #include "PscuHeader.h"
53 #include "TMultiGraph.h"
54 #include "physics/S4/S4Event.h"
55 #include "varDump/VarDumpEvent.h"
56 #include "varDump/VarDumpRecord.h"
57 #include "physics/trigger/TriggerEvent.h"
58
59 using namespace std;
60
61 void S4_QL(TString base, TString outDir, TString format, ULong_t DeltaT){ //DeltaT in minute
62
63 //------ Variables initialization ---------/
64 Int_t tmpSize;
65 ULong_t mintime, maxtime;
66 Int_t adcmax;
67 Int_t j=0;
68 Int_t S4_TRHadc;
69 Int_t S4_TRHmip;
70 TString filename = ((TObjString*)base.Tokenize('/')->Last())->GetString();
71 filename = ((TObjString*)filename.Tokenize('.')->First())->GetString();
72 char *trc;
73 TString trcs;
74 TString str;
75 TString trcsstot[31];
76 TString trcsstot2;
77 TString trcss[100]="";
78 TString trgconf[31]={"TOF1","TOF2","TOF3","TOF4","TOF5","TOF6","TOF7","S4","CALO","CALO-S4","TOF1-S4","TOF2-S4","TOF3-S4","TOF4-S4","TOF5-S4","TOF6-S4","TOF7-S4","TOF1-CALO","TOF2-CALO","TOF3-CALO","TOF4-CALO","TOF5-CALO","TOF6-CALO","TOF7-CALO","TOF1-CALO-S4","TOF2-CALO-S4","TOF3-CALO-S4","TOF4-CALO-S4","TOF5-CALO-S4","TOF6-CALO-S4","TOF7-CALO-S4",};
79 stringstream oss, oss1, oss2, oss3, s4soglia, buffer, conftrig, noentries;
80 ULong_t lastime, firstime, obt1;
81 Int_t vardumpentries = 0;
82 char S4_TRH[10];
83 char S4_TRH2[10];
84 Int_t trigconf = 0;
85 Long64_t nevents;
86 string title;
87 double obmin=0.;
88 double obmax=0.;
89 double obt;
90 double s4rate;
91 int nbin=0;
92 //------to open headerfile, s4file, Trigfile and vardumpfile---------------------------/
93 TFile *file =new TFile(base.Data()) ;
94 if (!file){
95 printf("No such file in the directory has been found");
96 return;
97 }
98
99 TTree *VarDumpTr = (TTree*)file->Get("VarDump");
100 TTree *PhysicsTr = (TTree*)file->Get("Physics");
101
102 TBranch *S4Br = PhysicsTr->GetBranch("S4");
103 TBranch *TriggerBr = PhysicsTr->GetBranch("Trigger");
104 TBranch *headBr = PhysicsTr->GetBranch("Header");
105 TBranch *VarDumpBr = VarDumpTr->GetBranch("VarDump");
106 TBranch *headVarDumpBr = VarDumpTr->GetBranch("Header");
107
108 pamela::S4::S4Event *s4e = 0;
109 pamela::EventHeader *eh = 0;
110 pamela::PscuHeader *ph = 0;
111 pamela::trigger::TriggerEvent *trige = 0;
112 pamela::VarDumpEvent *vde = 0;
113 pamela::VarDumpRecord *vdr = 0;
114
115 PhysicsTr->SetBranchAddress("S4", &s4e);
116 PhysicsTr->SetBranchAddress("Header", &eh);
117 PhysicsTr->SetBranchAddress("Trigger", &trige);
118 VarDumpTr->SetBranchAddress("VarDump", &vde);
119
120 nevents = headBr->GetEntries();
121 const Int_t sizetot = nevents;
122
123 //----------- If nevents < = 0 ---------------------------------/
124 if (nevents<=0) {
125 printf("nevents = %i \n", nevents);
126 printf(" \n");
127
128 TCanvas *canvas4 = new TCanvas("No entries", "No entries ", 400, 200);
129 canvas4->SetFillColor(10);
130 canvas4->cd();
131
132 TLatex *l = new TLatex();
133 l->SetTextAlign(12);
134 l->SetTextSize(0.15);
135 l->SetTextColor(2);
136 noentries.str("");
137 noentries<< "S4_QL:";
138 l->DrawLatex(0.05, 0.7, noentries.str().c_str());
139 noentries.str("");
140 noentries<< "No entries for this files";
141 l->DrawLatex(0.05, 0.5, noentries.str().c_str());
142
143 if (outDir == "./") {
144 oss.str("");
145 oss << filename.Data() << "_S4_QL." << format.Data();
146 } else {
147 oss.str("");
148 oss << outDir.Data() << filename.Data() << "_S4_QL." << format.Data();
149 }
150
151 canvas4->Update();
152 canvas4->SaveAs(oss.str().c_str());
153
154 return;
155 }
156
157 for (Int_t i = 0; i < nevents; i++){
158 headBr->GetEntry(i);
159 S4Br->GetEntry(i);
160 TriggerBr->GetEntry(i);
161 if (s4e->unpackError == 1 || s4e->S4_DATA==0) continue;
162 //-------to set max adc value------------------------------------------//
163 if (i==0){
164 adcmax=s4e->S4_DATA;
165 }
166 //-------to set trigger configuration---------------------------------//
167 if ((s4e->S4_DATA)> adcmax) adcmax=(s4e->S4_DATA);
168 trigconf = trige->trigconf;
169 trc = 0;
170 if ( trigconf & (1<<0) ) trc = "TOF1";
171 if ( trigconf & (1<<1) ) {
172 if (trc==0) trc= "TOF2";
173 else trc = Form("%s-TOF2",trc);
174 }
175 if ( trigconf & (1<<2) ) {
176 if (trc==0) trc= "TOF3";
177 else trc = Form("%s-TOF3",trc);
178 }
179 if ( trigconf & (1<<3) ) {
180 if (trc==0) trc= "TOF4";
181 else trc = Form("%s-TOF4",trc);
182 }
183 if ( trigconf & (1<<4) ) {
184 if (trc==0) trc= "TOF5";
185 else trc = Form("%s-TOF5",trc);
186 }
187 if ( trigconf & (1<<5) ) {
188 if (trc==0) trc= "TOF6";
189 else trc = Form("%s-TOF6",trc);
190 }
191 if ( trigconf & (1<<6) ) {
192 if (trc==0) trc= "TOF7";
193 else trc = Form("%s-TOF7",trc);
194 }
195 if ( trigconf & (1<<7) ) {
196 if (trc==0) trc= "S4";
197 else trc = Form("%s-S4",trc);
198 }
199 if ( trigconf & (1<<8) ) {
200 if (trc==0) trc= "CALO";
201 else trc = Form("%s-CALO",trc);
202 }
203 if ( trigconf & (1<<9) ) {
204 if (trc==0) trc= "CALIB_ON";
205 else trc = Form("%s-CALIB_ON",trc);
206 }
207 trcs = "";
208 trcs = trc;
209 if (i==0){
210 trcss[0]=trcs;
211 j=j+1;
212 }
213 if (trcs!=trcss[j-1] && i>0){
214 trcss[j]=trcs;
215 j=j+1;
216 }
217 }
218 for(Int_t k=0; k<31; k++){
219 for(Int_t p=0; p<j ; p++){
220 if(trgconf[k]==trcss[p]){
221 trcsstot[k]=trgconf[k];
222 }
223 }
224 }
225 for(Int_t k=0; k<31; k++){
226 if(trcsstot[k]!= "")
227 trcsstot2=trcsstot2+"/"+ trcsstot[k];
228 }
229 //---------to search S4 threshold and convert it into char (hexadecimal value)--------//
230 vardumpentries = VarDumpBr->GetEntries();
231 if(vardumpentries==0){
232 cout<<"TTree VarDump: Entries = 0"<<"\n"<<" S4 Threshold not defined"<<"\n";
233 str= " Not defined";
234 }
235 else{
236 for(Int_t k=0; k< vardumpentries; k++){
237 VarDumpBr->GetEntry(k);
238 vdr = (pamela::VarDumpRecord*)vde->Records->At(118);
239 S4_TRHadc=((vdr->VAR_VALUE)/256); //trh (number of ADC channels)
240 S4_TRHmip=S4_TRHadc*2;
241
242 //-----to convert decimal value to TString of S4 threshold value-------//
243 if(k==0){
244 sprintf(S4_TRH, "%d" , S4_TRHmip);
245 str= S4_TRH;
246 }else{
247 sprintf(S4_TRH2, "%d" , S4_TRHmip);
248 if(!strcmp(S4_TRH, S4_TRH2)) continue;
249 str= str+"/"+S4_TRH2;
250 }
251 sprintf(S4_TRH, "%d" , S4_TRHmip);
252 }
253 }
254
255 //********************** First Histogram ************************************//
256 TH1F *h1 = new TH1F("All events", "S4 distribution for file: "+filename, adcmax, 10, adcmax+(adcmax/10));
257 TH1F *h2 = new TH1F("S4 Triggered events", "S4 distribution for file: "+filename, adcmax, 10, adcmax+(adcmax/10));
258
259 //********************** Second Histogram ************************************//
260 for (Int_t i = 0; i < nevents; i++){
261 headBr->GetEntry(i);
262 ph = eh->GetPscuHeader();
263 obt1 = ph->GetOrbitalTime();
264 if(obt1 <= firstime) firstime=obt1;
265 if(obt1 >= lastime) lastime=obt1;
266 }
267 obmin=firstime;
268 obmax=lastime;
269 const ULong_t nint=(((lastime-firstime)/(DeltaT*60000)));
270 const Int_t size = (Int_t)(nint+1);
271
272 TH1F *Allev = new TH1F("Mean signal from S4-all events", filename+": Mean signal from S4 (all triggered events)", size, obmin, obmax);
273 TH1F *Alltime = new TH1F("Mean signal from S4-all events", filename+": Mean signal from S4 (all triggered events)", size, obmin, obmax);
274 TH1F *S4ev = new TH1F("Mean signal from S4-triggered events", filename+": Mean signal from S4 (only s4 triggered events)", size, obmin, obmax);
275 TH1F *S4time = new TH1F("Mean signal from S4-triggered events", filename+": Mean signal from S4 (only s4 triggered events)", size, obmin, obmax);
276
277 //********************** Third Histogram ************************************//
278 nbin=sizetot;
279 title="";
280 title=filename+": S4 rate from Trigger Packet";
281
282 const ULong_t nint3=(lastime-firstime);
283 const Int_t size3 = (Int_t)((nint3)/100);
284
285 TH1F *rate= new TH1F(title.c_str(), title.c_str(), nint3, obmin, obmax);
286 TH1F *rateline= new TH1F(filename+". S4 rate from Trigger Packet: mean value over 100 events", filename+". S4 rate from Trigger Packet: mean value over 100 events", size3, obmin, obmax);
287
288 //------------------------------------------------------------------------------------------------------------------
289 //------- fill histograms ---------//
290 Int_t n=0, p=0;
291 for (Int_t i = 0; i < nevents; i++){
292 TriggerBr->GetEntry(i);
293 headBr->GetEntry(i);
294 S4Br->GetEntry(i);
295 ph = eh->GetPscuHeader();
296 obt = ph->GetOrbitalTime();
297 if (s4e->unpackError == 1 && (s4e->S4_DATA)==0) continue;
298 s4rate= trige->s4calcount[0];
299 rate->Fill(obt, s4rate);
300 h1->Fill(s4e->S4_DATA);
301 Allev->Fill(obt,s4e->S4_DATA);
302 Alltime->Fill(obt);
303 if ((trige->patterntrig[0] == 0) && (trige->patterntrig[1] != 0) &&(trige->patterntrig[2] == 0) && (trige->patterntrig[3] == 0) && (trige->patterntrig[4] == 0) && (trige->patterntrig[5] == 0)){
304 h2->Fill(s4e->S4_DATA);
305 S4ev->Fill(obt,s4e->S4_DATA);
306 S4time->Fill(obt);
307 p=p+1;
308 }
309 }
310 Int_t kk=0;
311 while (kk < nevents){
312 obt=0;
313 s4rate=0;
314 for(Int_t jj=kk; jj< (kk+100); jj++){
315 TriggerBr->GetEntry(jj);
316 headBr->GetEntry(jj);
317 obt = obt+(ph->GetOrbitalTime());
318 s4rate= s4rate+(trige->s4calcount[0]);
319 }
320 rateline->Fill((obt/100),(s4rate/100));
321 kk=kk+100;
322 }
323
324 //****************************** Canvases *******************************//
325 //------------------- First Canvas --------------------------------//
326 TCanvas *canvas1 = new TCanvas("S4_QL_1", "S4 HISTO ", 1280, 1024);
327 canvas1->SetFillColor(10);
328 canvas1->Divide(1,2);
329
330 TPad *all= new TPad ("","", 0, 0, 1, 1);
331 all->SetFillColor(10);
332 TPad *s4 = new TPad ("s4","s4", 0, 1, 1, 0);
333 s4->SetFillColor(10);
334 TPad *pad = new TPad("pad","pad", .30,.80,.78,.90);
335 pad->SetFillColor(10);
336 TPad *pad2 = new TPad("pad2","pad2", .30,.80,.78,.90);
337 pad2->SetFillColor(10);
338 TLatex *l = new TLatex();
339 l->SetTextAlign(12);
340 l->SetTextSize(0.35);
341 l->SetTextColor(8);
342 s4soglia.str("");
343 s4soglia << "S4_THRESHOLD: "<<str.Data()<<" MIP";
344 TLatex *ll = new TLatex();
345 ll->SetTextAlign(12);
346 ll->SetTextSize(0.35);
347 conftrig.str("");
348 conftrig<<"TRIGGER CONF: "<<trcsstot2.Data();
349 if(adcmax<=100) adcmax=120;
350
351 canvas1->cd(1);
352 all->Draw();
353 all->cd();
354 h1->SetLineColor(kBlack);
355 h1->SetFillColor(kRed);
356 h1->GetXaxis()->SetTitle("ADC");
357 h1->GetXaxis()->CenterTitle();
358 h1->GetYaxis()->SetTitle("Number of Events");
359 h1->GetYaxis()->CenterTitle();
360 h1->Draw();
361 pad->Draw();
362 pad->cd();
363 l->DrawLatex(0.05, 0.65, s4soglia.str().c_str());
364 ll->DrawLatex(0.05, 0.15, conftrig.str().c_str());
365
366
367 canvas1->cd(2);
368 s4->Draw();
369 s4->cd();
370 h2->SetLineColor(kBlack);
371 h2->SetFillColor(2);
372 h2->GetXaxis()->SetTitle("ADC");
373 h2->GetXaxis()->CenterTitle();
374 h2->GetYaxis()->SetTitle("Number of Events");
375 h2->GetYaxis()->CenterTitle();
376 h2->Draw();
377 pad2->Draw();
378 pad2->cd();
379 l->DrawLatex(0.05, 0.65, s4soglia.str().c_str());
380 ll->DrawLatex(0.05, 0.15, conftrig.str().c_str());
381
382 if (h1->GetMaximum()>0){
383 all->SetLogy();
384 all->SetLogx();
385 }
386 if (h2->GetMaximum()>0){
387 s4->SetLogy();
388 s4->SetLogx();
389 }
390 canvas1->Update();
391
392
393 if (outDir == "./") {
394 oss1.str("");
395 oss1 << filename.Data() << "_S4_QL_1." << format.Data();
396 } else {
397 oss1.str("");
398 oss1 << outDir.Data() << filename.Data() << "_S4_QL_1." << format.Data();
399 }
400 canvas1->SaveAs(oss1.str().c_str());
401
402 //------------------- Second Canvas --------------------------------//
403 if (nint==0){
404 cout<<"Number of Time Interval = 0"<<"\n"<<" Set another value for input parameter DeltaTevtime "<<"\n"<<"RETURN"<<"\n";
405
406 TCanvas *canvas2 = new TCanvas("Number of time interval=0", "Number of time interval=0", 900, 300);
407 canvas2->SetFillColor(10);
408 canvas2->cd();
409
410 TLatex *l = new TLatex();
411 l->SetTextAlign(12);
412 l->SetTextSize(0.10);
413 l->SetTextColor(2);
414 noentries.str("");
415 noentries<< "S4_QL: Time evolution";
416 l->DrawLatex(0.05, 0.7, noentries.str().c_str());
417 noentries.str("");
418 noentries<< "Number of Time Interval = 0";
419 l->SetTextColor(1);
420 l->DrawLatex(0.05, 0.5, noentries.str().c_str());
421 noentries.str(" Set another value for input parameter DeltaTevtime");
422 noentries<< "";
423 l->DrawLatex(0.05, 0.3, noentries.str().c_str());
424
425 if (outDir == "./") {
426 oss.str("");
427 oss << filename.Data() << "_S4_QL_2." << format.Data();
428 } else {
429 oss.str("");
430 oss << outDir.Data() << filename.Data() << "_S4_QL_2." << format.Data();
431 }
432 canvas2->Update();
433 canvas2->SaveAs(oss.str().c_str());
434
435 return;
436 }
437 TCanvas *canvas2 = new TCanvas("S4_QL_2", "S4 - Time evolution", 1280, 1024);
438 canvas2->SetFillColor(10);
439 canvas2->Divide(1,2);
440
441 TPad *pad3 = new TPad("pad3","pad3", .7,.1,.88,.30);
442 pad3->SetFillColor(10);
443 TPad *pad4 = new TPad("pad4","pad4", .7,.1,.88,.30);
444 pad4->SetFillColor(10);
445
446 TLatex *l3 = new TLatex();
447 l3->SetTextAlign(12);
448 l3->SetTextSize(.3);
449 l3->SetTextColor(2);
450 buffer.str("");
451 buffer<<"DT : "<<DeltaT<<" min";
452
453 canvas2->cd(1);
454 Allev->SetStats(kFALSE);
455 Allev->SetMarkerColor(2);
456 Allev->SetMarkerSize(.5);
457 Allev->SetMarkerStyle(21);
458 Allev->GetXaxis()->SetTitle("Time ( OBT in ms )");
459 Allev->GetXaxis()->CenterTitle();
460 Allev->GetYaxis()->SetTitle("Mean signal ( ADC )");
461 Allev->GetYaxis()->CenterTitle();
462 Allev->SetMinimum(20);
463 Allev->Divide(Alltime);
464 Allev->Draw("p");
465 pad3->Draw();
466 pad3->cd();
467 l3->DrawLatex(0.05, 0.65, buffer.str().c_str());
468
469
470 canvas2->cd(2);
471 if(p != 0){
472 S4ev->SetStats(kFALSE);
473 S4ev->SetMarkerColor(2);
474 S4ev->SetMarkerSize(.5);
475 S4ev->SetMarkerStyle(21);
476 S4ev->GetXaxis()->SetTitle("Time ( OBT in ms )");
477 S4ev->GetXaxis()->CenterTitle();
478 S4ev->GetYaxis()->SetTitle("Mean signal ( ADC )");
479 S4ev->GetYaxis()->CenterTitle();
480 S4ev->SetMinimum(20);
481 S4ev->Divide(S4time);
482 S4ev->Draw("p");
483 pad4->Draw();
484 pad4->cd();
485 l3->DrawLatex(0.05, 0.65, buffer.str().c_str());
486 }else{
487 TLatex *l = new TLatex();
488 l->SetTextAlign(12);
489 l->SetTextSize(0.07);
490 l->SetTextColor(2);
491 noentries.str("");
492 noentries<< "S4_QL: Time evolution";
493 l->DrawLatex(0.06, 0.6, noentries.str().c_str());
494 noentries.str("");
495 noentries<< "No S4 triggered events.";
496 l->SetTextColor(1);
497 l->DrawLatex(0.06, 0.5, noentries.str().c_str());
498 }
499
500 if (outDir == "./") {
501 oss2.str("");
502 oss2 << filename.Data() << "_S4_QL_2." << format.Data();
503 } else {
504 oss2.str("");
505 oss2 << outDir.Data() << filename.Data() << "_S4_QL_2." << format.Data();
506 }
507 canvas2->SaveAs(oss2.str().c_str());
508
509 //------------------- Third Canvas --------------------------------//
510 TCanvas *canvas3 = new TCanvas("S4 - Rate", "S4 - Rate", 1280, 1024);
511 canvas3->SetFillColor(10);
512 canvas3->Divide(1,2);
513
514 canvas3->cd(1);
515 rate->SetStats(kFALSE);
516 rate->SetMarkerColor(2);
517 rate->SetMarkerSize(.5);
518 rate->SetMarkerStyle(21);
519 rate->GetXaxis()->SetTitle("Time ( OBT in ms )");
520 rate->GetXaxis()->CenterTitle();
521 rate->GetYaxis()->SetTitle("S4 rate (Hz)");
522 rate->GetYaxis()->CenterTitle();
523 if(rate->GetMaximum() > 1000) gPad->SetLogy();
524 rate->Draw("9p");
525
526 canvas3->cd(2);
527 rateline->SetStats(kFALSE);
528 rateline->SetMarkerColor(4);
529 rateline->SetMarkerSize(.5);
530 rateline->SetMarkerStyle(21);
531 rateline->GetXaxis()->SetTitle("Time ( OBT in ms )");
532 rateline->GetXaxis()->CenterTitle();
533 rateline->GetYaxis()->SetTitle("S4 rate (Hz)");
534 rateline->GetYaxis()->CenterTitle();
535 rateline->SetMaximum(rate->GetMaximum());
536 if(rateline->GetMaximum() > 1000) gPad->SetLogy();
537 rateline->Draw("9p");
538
539
540 if (outDir == "./") {
541 oss3.str("");
542 oss3 << filename.Data() << "_S4_QL_3." << format.Data();
543 } else {
544 oss3.str("");
545 oss3 << outDir.Data() << filename.Data() << "_S4_QL_3." << format.Data();
546 }
547 canvas3->SaveAs(oss3.str().c_str());
548 file->Close();
549 }
550
551
552 int main(int argc, char* argv[]){
553 TString path;
554 TString outDir ="./";
555 TString format ="jpg";
556 ULong_t DeltaT =1;
557
558
559 if (argc < 2){
560 printf("You have to insert at least the file to analyze \n");
561 printf("Try '--help' for more information. \n");
562 exit(1);
563 }
564
565 if (!strcmp(argv[1], "--help")){
566 printf( "Usage: ND_QL FILE [OPTION] \n");
567 printf( "\t --help Print this help and exit \n");
568 printf( "\t -outDir[path] Path where to put the output [default ./] \n");
569 printf( "\t -format[jpg|gif|ps] Format for output files [default 'jpg'] \n");
570 printf( "\t -DeltaT Time interval to control time evolution of acquired data in minute [default = 1 min] \n");
571 exit(1);
572 }
573
574 path=argv[1];
575
576 for (int i = 2; i < argc; i++){
577
578 if (!strcmp(argv[i], "-outDir")){
579 if (++i >= argc){
580 printf( "-outDir needs arguments. \n");
581 printf( "Try '--help' for more information. \n");
582 exit(1);
583 }
584 else{
585 outDir = argv[i];
586 continue;
587 }
588 }
589
590 if (!strcmp(argv[i], "-format")){
591 if (++i >= argc){
592 printf( "-format needs arguments. \n");
593 printf( "Try '--help' for more information. \n");
594 exit(1);
595 }
596 else{
597 format = argv[i];
598 continue;
599 }
600 }
601
602 if (!strcmp(argv[i], "-DeltaT")){
603 if (++i >= argc){
604 printf( "-DeltaT needs arguments. \n");
605 printf( "Try '--help' for more information. \n");
606 exit(1);
607 }
608 else{
609 DeltaT = atol(argv[i]);
610 continue;
611 }
612 }
613
614
615 }
616
617 S4_QL(argv[1], outDir, format, DeltaT);
618
619 }
620

  ViewVC Help
Powered by ViewVC 1.1.23