]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/AliFMDBaseDA.cxx
Add coloured boxes on DQM plots for easy reading.
[u/mrichter/AliRoot.git] / FMD / AliFMDBaseDA.cxx
1 /**************************************************************************
2  * Copyright(c) 2004, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /** @file    AliFMDBaseDA.cxx
17     @author  Hans Hjersing Dalsgaard <canute@nbi.dk>
18     @date    Wed Mar 26 11:30:45 2008
19     @brief   Base class for detector algorithms.
20 */
21 //
22 // This is the implementation of the (virtual) base class for the FMD
23 // detector algorithms(DA). It implements the creation of the relevant
24 // containers and handles the loop over the raw data. The derived
25 // classes can control the parameters and action to be taken making
26 // this the base class for the Pedestal, Gain and Physics DA.
27 //
28
29 #include "AliFMDBaseDA.h"
30 #include "AliRawReader.h"
31 #include "AliFMDDigit.h"
32 #include "AliFMDParameters.h"
33 #include "AliFMDRawReader.h"
34 #include "AliFMDCalibSampleRate.h"
35 #include "AliFMDCalibStripRange.h"
36 #include "AliLog.h"
37 #include "AliRawEventHeaderBase.h"
38 #include "AliFMDDigit.h"
39 #include <TClonesArray.h>
40 #include <TFile.h>
41 #include <TDatime.h>
42 #include <TSystem.h>
43 #include <TH2F.h>
44 #include <TStopwatch.h>
45 #include <TROOT.h>
46 #include <TPluginManager.h>
47 #include <iostream>
48 #include <fstream>
49
50 //_____________________________________________________________________
51 ClassImp(AliFMDBaseDA)
52 #if 0 
53 ; // Do not delete  - to let Emacs for mat the code
54 #endif
55
56 //_____________________________________________________________________
57 const char*
58 AliFMDBaseDA::GetStripPath(UShort_t det, 
59                            Char_t   ring, 
60                            UShort_t sec, 
61                            UShort_t str, 
62                            Bool_t   full) const
63 {
64   // Get the strip path 
65   // 
66   // Parameters 
67   //     det      Detector number
68   //     ring     Ring identifier 
69   //     sec      Sector number 
70   //     str      Strip number
71   //     full     If true, return full path 
72   // 
73   // Return 
74   //     The path
75   return Form("%s%sFMD%d%c[%02d,%03d]", 
76               (full ? GetSectorPath(det, ring, sec, full) : ""), 
77               (full ? "/" : ""), det, ring, sec, str);
78 }
79 //_____________________________________________________________________
80 const char*
81 AliFMDBaseDA::GetSectorPath(UShort_t det, 
82                             Char_t   ring, 
83                             UShort_t sec, 
84                             Bool_t   full) const
85 {
86   // Get the strip path 
87   // 
88   // Parameters 
89   //     det      Detector number
90   //     ring     Ring identifier 
91   //     sec      Sector number 
92   //     str      Strip number
93   //     full     If true, return full path 
94   // 
95   // Return 
96   //     The path
97   return Form("%s%sFMD%d%c[%02d]", 
98               (full ? GetRingPath(det, ring, full) : ""), 
99               (full ? "/" : ""), det, ring, sec);
100 }
101 //_____________________________________________________________________
102 const char*
103 AliFMDBaseDA::GetRingPath(UShort_t det, 
104                           Char_t   ring, 
105                           Bool_t   full) const
106 {
107   // Get the strip path 
108   // 
109   // Parameters 
110   //     det      Detector number
111   //     ring     Ring identifier 
112   //     sec      Sector number 
113   //     str      Strip number
114   //     full     If true, return full path 
115   // 
116   // Return 
117   //     The path
118   return Form("%s%sFMD%d%c", 
119               (full ? GetDetectorPath(det, full) : ""), 
120               (full ? "/" : ""), det, ring);
121 }
122 //_____________________________________________________________________
123 const char*
124 AliFMDBaseDA::GetDetectorPath(UShort_t det, 
125                               Bool_t   full) const
126 {
127   // Get the strip path 
128   // 
129   // Parameters 
130   //     det      Detector number
131   //     ring     Ring identifier 
132   //     sec      Sector number 
133   //     str      Strip number
134   //     full     If true, return full path 
135   // 
136   // Return 
137   //     The path
138   return Form("%s%sFMD%d", 
139               (full ? fDiagnosticsFilename.Data() : ""), 
140               (full ? ":/" : ""), det);
141 }
142
143 //_____________________________________________________________________
144 AliFMDBaseDA::AliFMDBaseDA() : 
145   TNamed(),
146   fDiagnosticsFilename("diagnosticsHistograms.root"),
147   fOutputFile(),
148   fConditionsFile(),
149   fSaveHistograms(kFALSE),
150   fMakeSummaries(kFALSE),
151   fDetectorArray(),
152   fPulseSize(10),
153   fPulseLength(10),
154   fRequiredEvents(0),
155   fCurrentEvent(0), 
156   fRunno(0),
157   fSummaries(0)
158 {
159   //Constructor
160   fSeenDetectors[0] = fSeenDetectors[1] = fSeenDetectors[2] = kFALSE;
161   fDetectorArray.SetOwner();
162   Rotate("conditions.csv", 3);
163   fConditionsFile.open("conditions.csv");
164 }
165 //_____________________________________________________________________
166 AliFMDBaseDA::AliFMDBaseDA(const AliFMDBaseDA & baseDA) : 
167   TNamed(baseDA),
168   fDiagnosticsFilename(baseDA.fDiagnosticsFilename.Data()),
169   fOutputFile(),
170   fConditionsFile(),
171   fSaveHistograms(baseDA.fSaveHistograms),
172   fMakeSummaries(baseDA.fMakeSummaries),
173   fDetectorArray(baseDA.fDetectorArray),
174   fPulseSize(baseDA.fPulseSize),
175   fPulseLength(baseDA.fPulseLength),
176   fRequiredEvents(baseDA.fRequiredEvents),
177   fCurrentEvent(baseDA.fCurrentEvent),
178   fRunno(baseDA.fRunno),
179   fSummaries(0)
180 {
181   //Copy constructor
182   fSeenDetectors[0] = baseDA.fSeenDetectors[0];
183   fSeenDetectors[1] = baseDA.fSeenDetectors[1];
184   fSeenDetectors[2] = baseDA.fSeenDetectors[2];
185
186   fDetectorArray.SetOwner();
187   
188 }
189
190
191 //_____________________________________________________________________
192 AliFMDBaseDA::~AliFMDBaseDA() 
193 {
194   //destructor
195 }
196
197 //_____________________________________________________________________
198 Bool_t AliFMDBaseDA::HaveEnough(Int_t nEvents) const
199 {
200   return nEvents > GetRequiredEvents();
201 }
202 //_____________________________________________________________________
203 UShort_t AliFMDBaseDA::GetProgress(Int_t nEvents) const
204 {
205   return UShort_t((nEvents *100)/ GetRequiredEvents());
206 }
207
208 //_____________________________________________________________________
209 void AliFMDBaseDA::Run(AliRawReader* reader) 
210 {
211   //Run the FMD DA
212   TFile* diagFile = 0;
213   // if (fSaveHistograms)
214   //  diagFile = TFile::Open(fDiagnosticsFilename.Data(),"RECREATE");
215   
216   reader->Reset();
217   fRunno = reader->GetRunNumber();
218
219   AliFMDRawReader* fmdReader  = new AliFMDRawReader(reader,0);
220   TClonesArray*    digitArray = new TClonesArray("AliFMDDigit",0);
221   
222   Bool_t sodread = kFALSE;
223   
224   for(Int_t i=0;i<3;i++) {
225     if (!reader->NextEvent()) { 
226       // Read Start-of-Run / Start-of-Files event
227       AliWarning(Form("Failed to read the %d%s event",
228                       i+1, (i == 0 ? "st" : (i == 1 ? "nd" : "rd"))));
229       break;
230     }
231     
232     UInt_t eventType = reader->GetType();
233     if(eventType == AliRawEventHeaderBase::kStartOfData || 
234        eventType == AliRawEventHeaderBase::kFormatError) { 
235       
236       WriteConditionsData(fmdReader);
237       Init();
238       sodread = kTRUE;
239       break;
240     }
241   }
242   
243   InitContainer(diagFile);
244   if (AliLog::GetDebugLevel("FMD","") >= 3) { 
245     fDetectorArray.ls();
246   }
247   
248   if(!sodread) 
249     AliWarning("No SOD event detected!");
250   
251   int lastProgress = 0;
252   
253   for(Int_t i = 0; i< 3;i++) fNEventsPerDetector[i] = 0;
254
255   for(Int_t n = 1; !HaveEnough(n); n++) {
256     if(!reader->NextEvent()) { n--; continue; }
257     SetCurrentEvent(n);
258     digitArray->Clear();
259     fmdReader->ReadAdcs(digitArray);
260     
261     Bool_t seen[] = { false, false, false };
262     for(Int_t i = 0; i<digitArray->GetEntriesFast();i++) {
263       AliFMDDigit* digit = static_cast<AliFMDDigit*>(digitArray->At(i));
264       UShort_t det = digit->Detector();
265       fSeenDetectors[det-1] = true;
266       seen[det-1]           = true;
267       FillChannels(digit);
268     }
269     
270     for(Int_t i = 0; i< 3;i++) 
271       if (seen[i]) (fNEventsPerDetector[i])++;
272       
273     FinishEvent();
274     
275     int progress = GetProgress(n);
276     if (progress <= lastProgress) continue;
277     lastProgress = progress;
278     std::cout << "Progress: " << lastProgress << " / 100 " << std::endl;
279     
280   }
281   
282   AliInfoF("Looped over %d events (%d,%d,%d)",GetCurrentEvent(),
283            fNEventsPerDetector[0], 
284            fNEventsPerDetector[1], 
285            fNEventsPerDetector[2]);
286   WriteHeaderToFile();
287   
288   for(UShort_t det=1;det<=3;det++) {
289     if (!fSeenDetectors[det-1]) continue;
290     std::cout << "FMD" << det << std::endl;
291     UShort_t firstRing = (det == 1 ? 1 : 0);
292     for (UShort_t ir = firstRing; ir < 2; ir++) {
293       Char_t   ring = (ir == 0 ? 'O' : 'I');
294       UShort_t nsec = (ir == 0 ? 40  : 20);
295       UShort_t nstr = (ir == 0 ? 256 : 512);
296
297       if (fMakeSummaries) MakeSummary(det, ring);
298
299       std::cout << " Ring " << ring << ": " << std::flush;
300       for(UShort_t sec =0; sec < nsec;  sec++)  {
301         for(UShort_t strip = 0; strip < nstr; strip++) {
302           Analyse(det,ring,sec,strip);
303         }
304         std::cout << '.' << std::flush;
305       }
306       // if(fSaveHistograms)
307       // diagFile->Flush();
308       std::cout << "done" << std::endl;
309     }
310   }
311   
312   if(fOutputFile.is_open()) {
313     fOutputFile.write("# EOF\n",6);
314     fOutputFile.close();
315   }
316   
317   Terminate(diagFile);
318     
319   if(fSaveHistograms ) {
320     diagFile = TFile::Open(fDiagnosticsFilename.Data(),"RECREATE");
321     fDetectorArray.Write("FMD", TObject::kSingleKey);
322     fSummaries.Write();
323     AliInfo("Closing diagnostics file - please wait ...");
324     // diagFile->Write();
325     diagFile->Close();
326     AliInfo("done");
327     
328   }
329 }
330 //_____________________________________________________________________
331
332 void AliFMDBaseDA::InitContainer(TDirectory* diagFile)
333 {
334   //Prepare container for diagnostics    
335   TDirectory* savDir   = gDirectory;
336
337   for(UShort_t det=1;det<=3;det++) {
338     TObjArray* detArray = new TObjArray(det == 1 ? 1 : 2);
339     detArray->SetOwner();
340     detArray->SetName(GetDetectorPath(det, false));
341     fDetectorArray.AddAtAndExpand(detArray,det);
342
343     TDirectory* detDir = 0;
344     if (diagFile) {
345       diagFile->cd();
346       detDir = diagFile->mkdir(GetDetectorPath(det, kFALSE));
347     }
348
349     UShort_t FirstRing = (det == 1 ? 1 : 0);
350     for (UShort_t ir = FirstRing; ir < 2; ir++) {
351       Char_t   ring = (ir == 0 ? 'O' : 'I');
352       UShort_t nsec = (ir == 0 ? 40  : 20);
353       UShort_t nstr = (ir == 0 ? 256 : 512);
354       TObjArray* ringArray = new TObjArray(nsec);
355       ringArray->SetOwner();
356       ringArray->SetName(GetRingPath(det, ring, false));
357       detArray->AddAtAndExpand(ringArray,ir);
358
359
360       TDirectory* ringDir = 0;
361       if (detDir) { 
362         detDir->cd();
363         ringDir = detDir->mkdir(GetRingPath(det,ring, kFALSE));
364       }
365       
366
367       for(UShort_t sec =0; sec < nsec;  sec++)  {
368         TObjArray* sectorArray = new TObjArray(nstr);
369         sectorArray->SetOwner();
370         sectorArray->SetName(GetSectorPath(det, ring, sec, false));
371         ringArray->AddAtAndExpand(sectorArray,sec);
372
373
374         TDirectory* secDir = 0;
375         if (ringDir) { 
376           ringDir->cd();
377           secDir = ringDir->mkdir(GetSectorPath(det, ring, sec, kFALSE));
378         }
379         
380         for(UShort_t strip = 0; strip < nstr; strip++) {
381           if (secDir) { 
382             secDir->cd();
383             secDir->mkdir(GetStripPath(det, ring, sec, strip, kFALSE));
384           }
385           TObjArray* stripArray = new TObjArray(0);
386           stripArray->SetOwner(true);
387           stripArray->SetName(GetStripPath(det, ring, sec, strip, false));
388           sectorArray->AddAtAndExpand(stripArray, strip);
389           AddChannelContainer(stripArray, det, ring, sec, strip);
390         }
391         AddSectorSummary(sectorArray, det, ring, sec, nstr);
392       }
393     }
394   }
395   savDir->cd();
396 }
397
398 //_____________________________________________________________________ 
399 void AliFMDBaseDA::WriteConditionsData(AliFMDRawReader* fmdReader) 
400 {
401   //Write the conditions data to file
402   AliFMDParameters* pars       = AliFMDParameters::Instance();
403   fConditionsFile.write(Form("# %s \n",pars->GetConditionsShuttleID()),14);
404   TDatime now;
405   fConditionsFile << "# This file created from run number " << fRunno 
406                   << " at " << now.AsString() << std::endl;
407   
408   AliFMDCalibSampleRate* sampleRate = new AliFMDCalibSampleRate();
409   AliFMDCalibStripRange* stripRange = new AliFMDCalibStripRange();
410   
411   fmdReader->ReadSODevent(sampleRate,stripRange,fPulseSize,fPulseLength,
412                           fSeenDetectors);
413
414   sampleRate->WriteToFile(fConditionsFile, fSeenDetectors);
415   stripRange->WriteToFile(fConditionsFile, fSeenDetectors);
416
417  
418   // Zero Suppresion
419   
420   // Strip Range
421   
422   fConditionsFile.write("# Gain Events \n",15);
423   
424   for(UShort_t det=1; det<=3;det++) {
425     if (!fSeenDetectors[det-1]) { 
426       continue;
427     }
428     UShort_t firstring = (det == 1 ? 1 : 0);
429     for(UShort_t iring = firstring; iring <=1;iring++) {
430       Char_t ring = (iring == 1 ? 'I' : 'O');
431       for(UShort_t board =0 ; board <=1; board++) {
432         
433         Int_t idx = GetHalfringIndex(det,ring,board);
434         
435         fConditionsFile << det                     << ','
436                         << ring                    << ','
437                         << board                   << ','
438                         << fPulseLength.At(idx)    << "\n";
439         
440       }
441     }
442   }
443   
444   fConditionsFile.write("# Gain Pulse \n",14);
445   
446   for(UShort_t det=1; det<=3;det++) {
447     if (!fSeenDetectors[det-1]) { 
448       continue;
449     }
450     UShort_t firstring = (det == 1 ? 1 : 0);
451     for(UShort_t iring = firstring; iring <=1;iring++) {
452       Char_t ring = (iring == 1 ? 'I' : 'O');
453       for(UShort_t board =0 ; board <=1; board++) {
454         
455         Int_t idx = GetHalfringIndex(det,ring,board);
456         
457         fConditionsFile << det                     << ','
458                         << ring                    << ','
459                         << board                   << ','
460                         << fPulseSize.At(idx)      << "\n";
461         
462       }
463     }
464   }
465   // sampleRate->WriteToFile(std::cout, fSeenDetectors);
466   // stripRange->WriteToFile(std::cout, fSeenDetectors);
467
468   if(fConditionsFile.is_open()) {
469     
470     fConditionsFile.write("# EOF\n",6);
471     fConditionsFile.close();
472     
473   }
474   
475 }
476 //_____________________________________________________________________ 
477 Int_t AliFMDBaseDA::GetHalfringIndex(UShort_t det, Char_t ring, 
478                                      UShort_t board) const 
479 {
480   // Get the index corresponding to a half-ring 
481   // 
482   // Parameters: 
483   //   det    Detector number 
484   //   ring   Ring identifier 
485   //   board  Board number 
486   //
487   // Return 
488   //   Internal index of the board 
489   UShort_t iring  =  (ring == 'I' ? 1 : 0);
490   
491   Int_t index = (((det-1) << 2) | (iring << 1) | (board << 0));
492   
493   return index-2;
494   
495 }
496 //_____________________________________________________________________ 
497 void AliFMDBaseDA::Rotate(const char* base, int max) const
498 {
499   // 
500   // Rotate a set of files.   base is the basic name of the files.
501   // If the file base.max exists it is removed. 
502   // If the file base.n exists (where n < max) it is renamed to
503   // base.(n-1).  
504   // If the file base exists, it is renamed to base.1 
505   //
506   // Parameters:
507   //   base Base name of the files
508   //   max  Maximum number to keep (minus one for the current).
509
510   // Note:  TSystem::AccessPathName returns false if the condition is
511   // fulfilled! 
512
513   // Check if we have base.max, and if so, remove it. 
514   TString testName(Form("%s.%d", base, max));
515   if (!gSystem->AccessPathName(testName.Data())) 
516     gSystem->Unlink(testName.Data());
517     
518   // Loop down from max-1 to 1 and move files 
519   for (int i = max-1; i >= 1; i--) { 
520     testName = Form("%s.%d", base, i);
521     if (!gSystem->AccessPathName(testName.Data())) {
522       TString newName(Form("%s.%d", base, i+1));
523       gSystem->Rename(testName.Data(), newName.Data());
524     }
525   }
526
527   // If we have the file base, rename it to base.1 
528   testName = Form("%s", base);
529   if (!gSystem->AccessPathName(testName.Data())){
530     TString newName(Form("%s.%d", base, 1));
531     gSystem->Rename(testName.Data(), newName.Data());
532   }
533 }
534
535 //_____________________________________________________________________ 
536 TH2*
537 AliFMDBaseDA::MakeSummaryHistogram(const char* prefix, const char* title, 
538                                    UShort_t d, Char_t r) 
539 {
540   // 
541   // Utility function for defining summary histograms 
542   // 
543   // Parameters:
544   //    det    Detector 
545   //    ring   Ring identifier 
546   //    prefix Histogram prefix 
547   //    title  Histogram title 
548   //
549   Int_t nX = ((d == 1 || r == 'I' || r == 'i') ?  20 :  40);
550   Int_t nY = ((d == 1 || r == 'I' || r == 'i') ? 512 : 256);
551   
552   TH2* ret = new TH2F(Form("%sFMD%d%c", prefix, d, r), 
553                       Form("%s for FMD%d%c", title, d, r), 
554                       nX, -0.5, nX-0.5, nY, -0.5, nY-0.5);
555   ret->SetXTitle("Sector #");
556   ret->SetYTitle("Strip #");
557   ret->SetDirectory(0);
558   // if (!fSummaries) fSummaries = new TObjArray;
559   fSummaries.Add(ret);
560   return ret;
561 }
562
563 //_____________________________________________________________________ 
564 TObjArray*
565 AliFMDBaseDA::GetDetectorArray(UShort_t det)
566 {
567   if (det < 1 || det > 3) { 
568     AliErrorF("Detector index %d out of bounds", det);
569     return 0;
570   }
571   return static_cast<TObjArray*>(fDetectorArray.At(det));
572 }
573 //_____________________________________________________________________ 
574 TObjArray*
575 AliFMDBaseDA::GetRingArray(UShort_t det, Char_t ring)
576 {
577   Int_t idx = (ring == 'O' || ring == 'o' ? 0 : 1);
578   TObjArray* array = GetDetectorArray(det);
579   if (!array) return 0;
580   array = static_cast<TObjArray*>(array->At(idx));
581   if (!array) AliErrorF("No ring array for FMD%d%c (%d)", det, ring, idx);
582   return array;
583 }
584 //_____________________________________________________________________ 
585 TObjArray*
586 AliFMDBaseDA::GetSectorArray(UShort_t det, Char_t ring, UShort_t sector)
587 {
588   TObjArray* array = GetRingArray(det,ring);
589   if (!array) return 0;
590   array = static_cast<TObjArray*>(array->At(sector));
591   if (!array) AliErrorF("No sector array for FMD%d%c[%02d]", det, ring, sector);
592   return array;
593 }
594 //_____________________________________________________________________ 
595 TObjArray*
596 AliFMDBaseDA::GetStripArray(UShort_t det, Char_t ring, 
597                             UShort_t sector, UShort_t strip)
598 {
599   TObjArray* array = GetSectorArray(det,ring,sector);
600   if (!array) return 0;
601   array = static_cast<TObjArray*>(array->At(strip));
602   if (!array) AliErrorF("No strip array for FMD%d%c[%02d,%03d]", 
603                         det, ring, sector, strip);
604   return array;
605 }
606
607 //=====================================================================
608 AliFMDBaseDA::Runner::Runner()
609   : fReader(0),
610     fDiagFile(""), 
611     fDiag(false)
612 {}
613
614 //_____________________________________________________________________ 
615 void
616 AliFMDBaseDA::Runner::AddHandlers()
617 {
618   gROOT->GetPluginManager()->AddHandler("TVirtualStreamerInfo",
619                                         "*",
620                                         "TStreamerInfo",
621                                         "RIO",
622                                         "TStreamerInfo()");
623   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", "Minuit", 
624                                         "TMinuitMinimizer",
625                                         "Minuit", 
626                                         "TMinuitMinimizer(const char *)");
627   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", 
628                                         "GSLMultiMin", 
629                                         "ROOT::Math::GSLMinimizer",
630                                         "MathMore", 
631                                         "GSLMinimizer(const char *)");
632   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", 
633                                         "GSLMultiFit", 
634                                         "ROOT::Math::GSLNLSMinimizer",
635                                         "MathMore", "GSLNLSMinimizer(int)");
636   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", 
637                                         "GSLSimAn", 
638                                         "ROOT::Math::GSLSimAnMinimizer",
639                                         "MathMore", 
640                                         "GSLSimAnMinimizer(int)");
641   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", 
642                                         "Linear", 
643                                         "TLinearMinimizer",
644                                         "Minuit", 
645                                         "TLinearMinimizer(const char *)");
646   gROOT->GetPluginManager()->AddHandler("ROOT::Math::Minimizer", 
647                                         "Fumili", 
648                                         "TFumiliMinimizer",
649                                         "Fumili", 
650                                         "TFumiliMinimizer(int)");
651 }
652 //_____________________________________________________________________ 
653 void
654 AliFMDBaseDA::Runner::ShowUsage(std::ostream& o, const char* progname)
655 {
656   o << "Usage: " << progname << " SOURCE [OPTIONS]\n\n"
657     << "Options:\n"
658     << "   -h,--help                Show this help\n"
659     << "   -d,--diagnostics[=FILE]  Write diagnostics to file\n"
660     << "   -D,--debug=LEVEL         Set the debug level\n\n"
661     << "SOURCE is one of\n"
662     << " * FILE.raw                Raw data file\n"
663     << " * FILE.root               ROOT'ified raw data file\n"
664     << " * collection://FILE.root  Event list in a ROOT file\n"
665     << " * collection://FILE       File containing list of ROOT files\n"
666     << " * ^FMD                    Monitor source\n"
667     << "There are other options too.  Check the AliRawReader docs\n"
668     << std::endl;
669 }
670  
671 //_____________________________________________________________________ 
672 namespace {
673   Bool_t ExtractValue(const TString& arg, TString& val)
674   {
675     val = "";
676     Int_t eq = arg.Index("=");
677     if (eq == kNPOS) return false;
678     
679     val = arg(eq+1, arg.Length()-eq-1);
680     return true;
681   }
682 }
683       
684 //_____________________________________________________________________ 
685 Int_t
686 AliFMDBaseDA::Runner::Init(int argc, char** argv)
687 {
688   AddHandlers();
689
690   // --- Process the command line ------------------------------------
691   TString source;
692   Int_t   debugLevel  = 0;
693   Bool_t  help        = false;
694   
695   for (int i = 1; i < argc; i++) { 
696     TString arg(argv[i]);
697     Bool_t  badOption   = false;
698     
699     if (arg[0] == '-') { // It's an option 
700       if (arg[1] == '-') { // It's a long option 
701         TString val;
702         if      (arg.EqualTo("--help"))     help = true; 
703         else if (arg.BeginsWith("--debug")) {
704           if (ExtractValue(arg, val))
705             debugLevel = val.Atoi();
706         }
707         else if (arg.BeginsWith("--diagnostics")) { 
708           fDiag = true;
709           if (ExtractValue(arg, val)) 
710             fDiagFile = val;
711         }
712         else badOption = true;
713       }
714       else { // Short option 
715         TString next(i < argc-1 ? argv[i+1] : "");
716         switch (arg[1]) { 
717         case 'h': help = true; break;
718         case 'd': fDiag = true; 
719           if (!next.IsNull() && next[0] != '-') {
720             fDiagFile = next;
721             i++;
722           }
723           break;
724         case 'D': 
725           if (!next.IsNull() && next[0] != '-') {
726             debugLevel = next.Atoi();
727             i++;
728           }
729           break;
730         default: badOption = true;
731         }
732       } // End of options
733       if (badOption) { 
734         std::cerr << argv[0] << ": Unknown option " << argv[i] 
735                   << std::endl;
736         return -1;
737       }
738     }
739     else { // source or compatibility debug level 
740       if (source.IsNull()) source = arg;
741       else                 debugLevel = arg.Atoi();
742     }
743   }
744   
745   // --- Check if help was requested ---------------------------------
746   if (help) { 
747     ShowUsage(std::cout, argv[0]);
748     return 1;
749   }
750
751   // --- Check if we have a source -----------------------------------
752   if (source.IsNull()) { 
753     std::cerr << "No source given" << std::endl;
754     return -2;
755   }
756
757   // --- Initialize various things -----------------------------------
758   AliFMDParameters::Instance()->Init(kFALSE,0);
759
760   //This will only work for FDR 1 data. When newer data becomes
761   //available the ! must be removed!
762   Bool_t old = kTRUE;
763   AliFMDParameters::Instance()->UseCompleteHeader(old);
764   
765   AliLog::EnableDebug(debugLevel > 0);
766   AliLog::SetModuleDebugLevel("FMD", debugLevel);
767
768   // --- Make our reader ---------------------------------------------
769   fReader = AliRawReader::Create(source);
770   if (!fReader) { 
771     std::cerr << "Failed to make raw reader for source \"" << source 
772               << "\"" << std::endl;
773     return -3;
774   }
775   return 0;
776 }
777
778 //_____________________________________________________________________ 
779 Int_t
780 AliFMDBaseDA::Runner::RunNumber() const
781
782   if (!fReader) return -1;
783   return fReader->GetRunNumber(); 
784 }
785
786 //_____________________________________________________________________ 
787 void
788 AliFMDBaseDA::Runner::Exec(AliFMDBaseDA& da)
789 {
790   TStopwatch timer;
791   timer.Start();
792
793   da.SetSaveDiagnostics(fDiag || !fDiagFile.IsNull());
794   if (!fDiagFile.IsNull()) da.SetDiagnosticsFilename(fDiagFile);
795
796   da.Run(fReader);
797
798   timer.Stop();
799   timer.Print();
800
801 #ifdef ALI_AMORE
802   try { 
803     amore::da::AmoreDA myAmore(amore::da::AmoreDA::kSender);
804
805     for (UShort_t det = 1; det <= 3; det++) {
806       if (!da.HasSeenDetector(det)) continue;
807       TObject* runNo = new TObject;
808       runNo->SetUniqueID(fReader->GetRunNumber());
809       myAmore.Send(Form("gainRunNoFMD%d", det), runNo);
810     }
811     
812     TIter     next(&da.GetSummaries());
813     TObject*  obj = 0;
814     while ((obj = next())) 
815       myAmore.Send(obj->GetName(), obj);
816     
817   }
818   catch (exception& e) {
819     cerr << "Failed to make AMORE instance: " << e.what() << endl;
820   }
821                                
822 #endif
823 }
824
825
826   
827   
828   
829 //_____________________________________________________________________ 
830
831 //_____________________________________________________________________ 
832 //
833 // EOF
834 //