]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/qa/QAPlotter.C
Merge branch 'master' of https://git.cern.ch/reps/AliRoot
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / qa / QAPlotter.C
1 /**
2  * @file   QAPlotter.C
3  * @author Christian Holm Christensen <cholm@nbi.dk>
4  * @date   Thu Nov 17 12:08:04 2011
5  * 
6  * @brief  Class to plot QA trends
7  * 
8  * @ingroup pwglf_forward_qa_scripts
9  * 
10  */
11 #ifndef __CINT__
12 # include "QABase.h"
13 # include <TTree.h>
14 # include <TGraph.h>
15 # include <TGraphErrors.h>
16 # include <TGraphAsymmErrors.h>
17 # include <TMultiGraph.h>
18 # include <TH1.h>
19 # include <TLegend.h>
20 # include <TCanvas.h>
21 # include <TLine.h>
22 # include <TArrow.h>
23 # include <TArrayI.h>
24 # include <TMath.h>
25 # include <TChain.h>
26 #else 
27 class QABase;
28 class QARing;
29 class TTree;
30 class TGraphAsymmErrors;
31 class TGraph;
32 class TGraphErrors;
33 class TMultiGraph;
34 class RingQuantity;
35 class TArrayI;
36 class TH1;
37 #endif
38
39
40 /**
41  * Class to plot QA trends 
42  * 
43  * @ingroup pwglf_forward_qa_scripts
44  */
45 struct QAPlotter : public QABase
46 {
47   /**
48    * Ring class 
49    */
50   struct Ring : public QARing
51   {
52     /** 
53      * Constuctor
54      * 
55      * @param d      Detector 
56      * @param r      Ring 
57      * @param useVar Use variance for errors (not min/max)
58      */
59     Ring(UShort_t d, Char_t r, Bool_t useVar=false)
60       : QARing(d, r),
61         fGChi2(0),
62         fGC(0),
63         fGDelta(0),
64         fGXi(0),
65         fGSigma(0),
66         fGLow(0),
67         fGSingles(0),
68         fGLoss(0),
69         fGBeta(0),
70         fGOccupancy(0),
71         fUseVar(useVar)
72     {
73       fGChi2           = new TGraphAsymmErrors;
74       fGC              = new TGraphAsymmErrors;
75       fGDelta          = new TGraphAsymmErrors;
76       fGXi             = new TGraphAsymmErrors;
77       fGSigma          = new TGraphAsymmErrors;
78       fGLow            = new TGraph;
79       fGSingles        = new TGraph;
80       fGLoss           = new TGraph;
81       fGBeta           = new TGraph;
82       fGOccupancy      = new TGraphAsymmErrors;
83
84       SetAtt(fGChi2,            "chi2",     "#LT#chi^{2}/#nu#GT");
85       SetAtt(fGC,               "c",        "Constant");
86       SetAtt(fGDelta,           "delta",    "#LT#Delta_{p}#GT");
87       SetAtt(fGXi,              "xi",       "#LT#xi#GT");
88       SetAtt(fGSigma,           "sigma",    "#LT#sigma#GT");
89       SetAtt(fGLow,             "low",      "# of low statistics bins");
90       SetAtt(fGSingles,         "singles",  "Fraction of single hits");
91       SetAtt(fGLoss,            "loss",     "Data lossed due to cuts");
92       SetAtt(fGOccupancy,       "occupancy","#LTOccupancy#GT");
93       SetAtt(fGBeta,            "beta",     "Correlation of methods");
94
95     }
96     /** 
97      * Static member function to get the ring color
98      * 
99      * @param d Detector number
100      * @param r Ring identifer 
101      * 
102      * @return Color
103      */
104     static Color_t RingColor(UShort_t d, Char_t r)
105     { 
106       return ((d == 1 ? kRed : (d == 2 ? kGreen : kBlue))
107               + ((r == 'I' || r == 'i') ? 2 : -3));
108     }
109     /** 
110      * Set graph attributes
111      * 
112      * @param g      Graph
113      * @param name   Name of graph
114      */
115     void SetAtt(TGraph* g, const char* name, const char* /*title=""*/) 
116     {
117       Color_t c = RingColor(fD, fR);
118       g->SetName(Form("FMD%d%c_%s", fD, fR, name));
119       // g->SetTitle(Form("FMD%d%c %s", fD, fR, 
120       // !title || title[0] == '\0' ? name : title));
121       Int_t marker = 20+(fD-1) + (fR == 'I' ? 0 : 4);
122       g->SetTitle(Form("FMD%d%c", fD, fR));
123       g->SetLineColor(c);
124       g->SetFillColor(c);
125       g->SetMarkerColor(c);
126       g->SetMarkerStyle(marker);
127       g->SetLineWidth(2);
128       switch (marker) { 
129       case 20: g->SetMarkerSize(1.2); break;
130       case 21: g->SetMarkerSize(1.2); break;
131       case 22: g->SetMarkerSize(1.3); break;
132       case 26: g->SetMarkerSize(1.1); break;
133       }
134     }
135     /** 
136      * Update a graph from a RingQuantity 
137      * 
138      * @param g      Graph to update
139      * @param q      Quantity 
140      * @param runNo  Run number 
141      */
142     void UpdateGraph(TGraphAsymmErrors* g, RingQuantity& q, 
143                      UInt_t runNo, UInt_t /* n */) 
144     {
145       Double_t y  = q.mean;
146       Double_t el = y-q.min;
147       Double_t eh = q.max-y;
148       if (fUseVar) { 
149         //Info("UpdateGraph", "Setting errors on %s to variance",g->GetName());
150         el = q.var; 
151         eh = q.var; 
152       }
153       if (TMath::Abs(y) < 1e-6) return;
154       Int_t    i  = g->GetN();
155       g->SetPoint(i, runNo, y);
156       g->SetPointError(i, 0, 0, el, eh);
157     }
158     /** 
159      * Update all graphs
160      * 
161      * @param n     Entry number (not used)
162      * @param runNo Run number 
163      * 
164      * @return true on success
165      */
166     Bool_t Update(UInt_t n, UInt_t runNo)
167     {
168       UpdateGraph(fGChi2,      *fChi2,          runNo, n);
169       UpdateGraph(fGC,         *fC,             runNo, n);
170       UpdateGraph(fGDelta,     *fDelta,         runNo, n);
171       UpdateGraph(fGXi,        *fXi,            runNo, n);
172       UpdateGraph(fGSigma,     *fSigma,         runNo, n);
173       UpdateGraph(fGOccupancy, *fOccupancy,     runNo, n);
174
175       fGLow->SetPoint(n, runNo, fFitStatus->nLow);
176       
177       if (fMerge->one > 1e-6)
178         fGSingles->SetPoint(fGSingles->GetN(), runNo, fMerge->one);
179       if (fCorrelation->beta > 1e-6)
180         fGBeta->SetPoint(fGBeta->GetN(), runNo, fCorrelation->beta);
181       if (-fDataLoss->full > 1e-6)
182         fGLoss->SetPoint(fGLoss->GetN(), runNo, -fDataLoss->full);
183       return true;
184     }
185     TGraphAsymmErrors* fGChi2;     // Graph of ELoss reduced chi-square
186     TGraphAsymmErrors* fGC;        // Graph of ELoss constant 
187     TGraphAsymmErrors* fGDelta;    // Graph of ELoss MPV
188     TGraphAsymmErrors* fGXi;       // Graph of ELoss Landau width
189     TGraphAsymmErrors* fGSigma;    // Graph of ELoss Gaus width
190     TGraph*            fGLow;      // Graph of bins with low statistics
191     TGraph*            fGSingles;  // Graph of fraction of singles
192     TGraph*            fGLoss;     // Graph of 'lost' data 
193     TGraph*            fGBeta;     // Graph of Poisson vs ELoss correlation
194     TGraphAsymmErrors* fGOccupancy;// Graph of mean occupancy              
195     Bool_t             fUseVar;    // Use variance 
196   };
197   /** 
198    * Constructor 
199    */
200   QAPlotter(const TString& dataType, 
201             Int_t          year, 
202             const TString& period, 
203             const TString& pass, 
204             Bool_t         useVar=true) 
205     : QABase(dataType, year, period, pass),
206       fNAccepted(0),
207       fVz(0),
208       fUseVar(useVar)
209   {
210     Info("QAPlotter", "Do we use variance? %s", fUseVar ? "yes" : "no");
211     fFMD1i = new Ring(1, 'I', useVar); 
212     fFMD2i = new Ring(2, 'I', useVar); 
213     fFMD2o = new Ring(2, 'O', useVar); 
214     fFMD3i = new Ring(3, 'I', useVar); 
215     fFMD3o = new Ring(3, 'O', useVar); 
216     fNAccepted = new TGraph;
217     fNAccepted->SetName("nAccepted");
218     fNAccepted->SetMarkerStyle(20);
219     fNAccepted->SetLineWidth(2);
220
221     fVz = new TGraphErrors;
222     fVz->SetName("vz");
223     fVz->SetMarkerStyle(20);
224     fVz->SetLineWidth(2);
225   }
226   /** 
227    * Add a file to be processed
228    * 
229    * @param filename Name of file 
230    */
231   void AddFile(const char* filename)
232   {
233     fFiles.Add(new TObjString(filename));
234   }
235   const char* GetUserInfo(TList* l, const char* name, const char* def) const
236   {
237     if (!l) {
238       Warning("GetUserInfo", "No user information list");
239       return def;
240     }
241     
242     TObject* o = l->FindObject(name);
243     if (!o) {
244       Warning("GetUserInfo", "User information %s not found", name);
245       l->ls();
246       return def;
247     }
248
249     Info("GetUserInfo", "Got user information %s=%s", name, o->GetTitle());
250     return o->GetTitle();
251   }
252   /** 
253    * Make a tree 
254    * 
255    * @param read If true, read from file 
256    * 
257    * @return True on success 
258    */
259   Bool_t MakeTree(bool read)
260   {
261     if (fFiles.GetEntriesFast() <= 0) return QABase::MakeTree(read);
262     if (fFiles.GetEntriesFast() == 1 && read) {
263       TFile* file = TFile::Open(fFiles.At(0)->GetName(), "READ");
264       if (file) {
265         fTree       = static_cast<TTree*>(file->Get("T"));
266         return true;
267       }
268     }
269     TChain* chain = new TChain("T", "T");
270     if (!chain->AddFileInfoList(&fFiles)) return false;
271     
272     fTree   = chain;
273     
274     return true;
275   }
276
277   /** 
278    * Run the job
279    * 
280    */
281   void Run()
282   {
283     Init(true);
284     if (!fTree) {
285       Error("Run", "No input tree");
286       return;
287     }
288     fFirst = 0xFFFFFFFF;
289     fLast  = 0;
290
291     UInt_t nEntries = fTree->GetEntries();
292     UInt_t j = 0;
293     fRuns.Set(nEntries);
294     Info("Run", "Got %d runs", nEntries);
295     TList* l = fTree->GetUserInfo();
296     fPeriod   = GetUserInfo(l, "period", "?");
297     fPass     = GetUserInfo(l, "pass",   "?");
298     fDataType = GetUserInfo(l, "type",   "?");
299       
300
301     for (UInt_t i = 0; i < nEntries; i++) {
302       fTree->GetEntry(i);
303
304       UInt_t run = fGlobal->runNo;
305       UInt_t nev = fGlobal->nAccepted;
306
307       fFirst = TMath::Min(run, fFirst);
308       fLast  = TMath::Max(run, fLast);
309       fRuns[i] = run;
310
311       Info("Run", "Got run %d with %d accepted events", run, nev);
312       fNAccepted->SetPoint(i, run, nev);
313       fVz->SetPoint(i, run, fGlobal->meanVz);
314       fVz->SetPointError(i, 0, fGlobal->sigmaVz);
315
316       if (nev <= 100) continue;
317       static_cast<Ring*>(fFMD1i)->Update(j, run);
318       static_cast<Ring*>(fFMD2i)->Update(j, run);
319       static_cast<Ring*>(fFMD2o)->Update(j, run);
320       static_cast<Ring*>(fFMD3i)->Update(j, run);
321       static_cast<Ring*>(fFMD3o)->Update(j, run);
322       j++;
323     }
324     
325     Plot();
326   }
327   /** 
328    * Plot results
329    * 
330    */
331   void Plot()
332   {
333     // fTeXName = Form("trend_%09d_%09d", fFirst, fLast);
334     TString title;
335     if (!fPeriod.IsNull() && !fPass.IsNull()) 
336       title.Form("QA trends for %s/%s runs %d --- %d", 
337                  fPeriod.Data(), fPass.Data(), fFirst, fLast);
338     else
339       title.Form("QA trends for runs %d --- %d", fFirst, fLast);
340     MakeCanvas(title);
341
342     CanvasTitle("# of accepted events");
343     fNAccepted->Draw("apl");
344     PutCanvasTitle("# of accepted events");
345     AddRuns(fNAccepted->GetHistogram(), "# of accepted events");
346
347     TLine* l = new TLine(fFirst, 100, fLast, 100);
348     l->SetLineColor(kRed+2);
349     l->SetLineStyle(2);
350     l->Draw();
351     PrintCanvas("nAccepted");
352
353     CanvasTitle("#LTv_{z}#GT");
354     fVz->Draw("apl");
355     PutCanvasTitle("Mean z coordinate of interaction point");
356     AddRuns(fVz->GetHistogram(), "#LTv_{z}#GT");
357     PrintCanvas("vz");
358
359     TMultiGraph* chi2           = new TMultiGraph;              
360     TMultiGraph* c              = new TMultiGraph;              
361     TMultiGraph* delta          = new TMultiGraph;              
362     TMultiGraph* xi             = new TMultiGraph;              
363     TMultiGraph* sigma          = new TMultiGraph;              
364     TMultiGraph* low            = new TMultiGraph;              
365     TMultiGraph* singles        = new TMultiGraph;      
366     TMultiGraph* loss           = new TMultiGraph;
367     TMultiGraph* occ    = new TMultiGraph;      
368     TMultiGraph* beta           = new TMultiGraph;              
369     chi2        ->SetName("chi2");             
370     c           ->SetName("c");                
371     delta       ->SetName("delta");            
372     xi          ->SetName("xi");               
373     sigma       ->SetName("sigma");            
374     low         ->SetName("low");              
375     singles     ->SetName("singles");    
376     loss        ->SetName("loss");
377     beta        ->SetName("beta");          
378     occ ->SetName("occupancy");
379
380
381     AddToMulti(fFMD1i,chi2, c, delta, xi, sigma, low, singles, loss, beta, occ);
382     AddToMulti(fFMD2i,chi2, c, delta, xi, sigma, low, singles, loss, beta, occ);
383     AddToMulti(fFMD2o,chi2, c, delta, xi, sigma, low, singles, loss, beta, occ);
384     AddToMulti(fFMD3i,chi2, c, delta, xi, sigma, low, singles, loss, beta, occ);
385     AddToMulti(fFMD3o,chi2, c, delta, xi, sigma, low, singles, loss, beta, occ);
386
387     PlotMulti(chi2,     "#LT#chi^{2}/#nu#GT from #Delta fits", true);
388     PlotMulti(c,        "#LTc#GT from #Delta fits");
389     PlotMulti(delta,    "#LT#Delta_{p}#GT from #Delta fits");
390     PlotMulti(xi,       "#LT#xi#GT from #Delta fits");
391     PlotMulti(sigma,    "#LT#sigma#GT from #Delta fits");
392     PlotMulti(low,      "Bins with too low statistics");
393     PlotMulti(singles,  "Fraction of single hits");
394     PlotMulti(loss,     "% of hits 'lost' due to merging+cuts");
395     PlotMulti(occ,      "#LTOccupancy#GT [%]", true);
396     PlotMulti(beta,     "Correlation of methods");
397
398     fStore->cd();
399     fNAccepted->Write();
400     chi2->Write();
401     c->Write();
402     delta->Write();
403     sigma->Write();
404     low->Write();
405     singles->Write();
406     loss->Write();
407     occ->Write();
408     beta->Write();
409
410     Close(false); // Do not delete PNGs
411   }
412   /** 
413    * Add graphs from a ring to the multi graphs
414    * 
415    * @param qr        Ring
416    * @param chi2      ELoss reduced chi-square     
417    * @param c         ELoss constant               
418    * @param delta     ELoss MPV            
419    * @param xi        ELoss Landau width           
420    * @param sigma     ELoss Gaus width             
421    * @param low       bins with low statistics     
422    * @param singles   fraction of singles          
423    * @param loss      'lost' data                  
424    * @param beta      Poisson vs ELoss correlation
425    * @param occupancy mean occupancy              
426    */
427   void AddToMulti(QARing* qr, 
428                   TMultiGraph* chi2,
429                   TMultiGraph* c,
430                   TMultiGraph* delta,
431                   TMultiGraph* xi,
432                   TMultiGraph* sigma,
433                   TMultiGraph* low,
434                   TMultiGraph* singles,
435                   TMultiGraph* loss,
436                   TMultiGraph* beta,
437                   TMultiGraph* occupancy)
438   {
439     Ring* r = static_cast<Ring*>(qr);
440     chi2        ->Add(r->fGChi2);
441     c           ->Add(r->fGC);
442     delta       ->Add(r->fGDelta);
443     xi          ->Add(r->fGXi);
444     sigma       ->Add(r->fGSigma);
445     low         ->Add(r->fGLow);
446     singles     ->Add(r->fGSingles);
447     loss        ->Add(r->fGLoss);
448     occupancy   ->Add(r->fGOccupancy);
449     beta        ->Add(r->fGBeta);
450   }
451   /** 
452    * Plot a multi-graph
453    * 
454    * @param mg     Multi graph
455    * @param title  Title
456    * @param logy   If true, make @f$\log@f$ scale 
457    */
458   void PlotMulti(TMultiGraph* mg, const char* title, Bool_t logy=false)
459   {
460     CanvasTitle(title);
461     // fCanvas->SetBottomMargin(.15);
462     fCanvas->SetLeftMargin(.08);
463     fCanvas->SetTopMargin(.1);
464     fCanvas->SetLogy(logy);
465     // fCanvas->SetRightMargin(.2);
466     mg->Draw("apl");
467     
468     TH1*     h   = mg->GetHistogram();
469     Double_t max = h->GetMaximum();
470     Double_t min = h->GetMinimum();
471     if (h->GetMinimum() == 0) {
472       min = min - .1*(max-min);
473       h->SetMinimum(min);
474     }
475     Int_t    x1  = h->GetXaxis()->GetXmin();
476     Int_t    x2  = h->GetXaxis()->GetXmax();
477     
478     TLegend* l = new TLegend(.1, .91, .97, .95);
479     l->SetNColumns(5);
480     l->SetFillColor(0);
481     l->SetFillStyle(0);
482     l->SetBorderSize(0);
483     l->SetTextFont(42);
484     TIter next(mg->GetListOfGraphs());
485     mg->GetListOfGraphs();
486     TGraph* g = 0;
487
488     // Get the runs we have here 
489     TArrayI runs(fRuns.GetSize());
490     runs.Reset(INT_MAX);
491     Int_t   j = 0;
492     while ((g = static_cast<TGraph*>(next()))) { 
493       l->AddEntry(g, g->GetTitle(), "lp");
494       Double_t* xs = g->GetX();
495       Int_t     n  = g->GetN();
496       for (Int_t i = 0; i < n; i++) {
497         if (FindRun(runs, Int_t(xs[i])) >= 0) continue;
498         runs.SetAt(xs[i], j++);
499       }
500       Double_t ymean = g->GetMean(2);
501       Double_t xh    = x2 - .03 * Double_t(x2-x1); 
502       TLine* lm = new TLine(x1, ymean, xh, ymean);
503       lm->SetLineColor(g->GetLineColor());
504       lm->SetLineStyle(2);
505       lm->SetLineWidth(1);
506       lm->Draw();
507       TLatex* al = new TLatex(xh, ymean, g->GetTitle());
508       al->SetTextColor(g->GetLineColor());
509       al->SetTextFont(42);
510       al->SetTextSize(.02);
511       al->SetTextAlign(12);
512       al->Draw();
513     }
514     l->Draw();
515
516     TList areas;
517     areas.SetOwner();
518     AddRuns(h, title, &runs, &areas);
519
520     PrintCanvas(mg->GetName(), &areas);
521   }
522   /** 
523    * Find a run 
524    * 
525    * @param runs List of runs 
526    * @param run  Run to find
527    * 
528    * @return Index of run in run list
529    */
530   Int_t FindRun(const TArrayI& runs, Int_t run) 
531   {
532     std::sort(&(runs.fArray[0]), &(runs.fArray[runs.GetSize()]));
533     Int_t idx = TMath::BinarySearch(runs.GetSize(), runs.fArray, run);
534     if (idx >= runs.GetSize() || idx < 0 || runs[idx] != run) return -1;
535     return idx;
536   }
537   /** 
538    * Add run labels at appropriate places on the plot
539    * 
540    * @param h      Frame histogram
541    * @param title  Title 
542    * @param runs   List of runs, if any
543    * @param areas  Other areas 
544    */
545   void AddRuns(TH1* h, const char* title, TArrayI* runs=0, 
546                TList* areas=0)
547   {
548     h->GetXaxis()->SetNoExponent();
549     // h->GetXaxis()->SetTitleOffset(1);
550     TString ytitle(title);
551     if (fUseVar) ytitle.Append(" (errors: variance)");
552     else         ytitle.Append(" (errors: min/max)");
553     h->SetYTitle(ytitle.Data());
554     h->SetXTitle("Run #");
555     // Info("AddRuns", "%s: %s vs %s", h->GetName(), 
556     //      h->GetXaxis()->GetTitle(), h->GetYaxis()->GetTitle());
557   
558     Int_t    r1  = h->GetXaxis()->GetXmin();
559     Int_t    r2  = h->GetXaxis()->GetXmax();
560     Double_t lx  = 0;
561     Double_t tx  = .025; // (r2 - r1) / 18;
562     Double_t wx  = 1 - fCanvas->GetLeftMargin() - fCanvas->GetRightMargin();
563     Double_t dy  = .025;
564     Double_t y   = fCanvas->GetBottomMargin()+dy;
565     for (Int_t i = 0; i < fRuns.GetSize(); i++) {
566       Int_t    r = fRuns[i];
567       Double_t x = fCanvas->GetLeftMargin() + wx*Double_t(r-r1)/(r2-r1);
568
569       // Skip runs out of range 
570       if (r < r1 || r > r2) continue;
571
572       // Skip runs not in the graphs 
573       if (runs) {
574         if (FindRun(*runs, r) < 0) { 
575           continue;
576         }
577       }
578           
579
580       if (TMath::Abs(x - lx) < tx) y += dy;
581       else                         y =  fCanvas->GetBottomMargin() + dy;
582
583
584       // Info(h->GetName(), "%6d (x,y)=(%12f,%12f) |lx-x|=|%f-x|=%f", 
585       //      r, x, y, lx, TMath::Abs(lx-x));
586       lx = x;
587       
588       Double_t* xa = fNAccepted->GetX();
589       Int_t     na = fNAccepted->GetN();
590       Int_t idx = TMath::BinarySearch(na, xa, Double_t(r));
591       Color_t color = kBlue+3;
592       if (idx >= 0  && idx < na  && r == xa[idx] && 
593           fNAccepted->GetY()[idx] < 10000) 
594         color = kRed+3;
595
596       TLatex* ll = new TLatex(x, y, Form("%d", r));
597       ll->SetNDC();
598       ll->SetTextAlign(21);
599       ll->SetTextSize(0.02);
600       ll->SetTextColor(color);
601       ll->SetTextFont(42);
602       // ll->SetTextAngle(90);
603       ll->Draw();
604       TLine* tl = new TLine(x, y, x, 1-fCanvas->GetTopMargin());
605       tl->SetBit(TLine::kLineNDC);
606       tl->SetLineStyle(3);
607       tl->SetLineColor(color);
608       tl->Draw();
609
610       if (!areas) continue;
611       
612       TObjString* area = new TObjString;
613       TString&    spec = area->String();
614       spec.Form("<span style=\"left: %d%%; bottom: %d%%;\" "
615                 "onClick='window.location=\"%09d/index.html\"' "
616                 "onMouseOver='this.style.cursor=\"pointer\"' "
617                 "alt=\"%d\""
618                 ">%d</span>",
619                 UInt_t(100*x)-2, UInt_t(100*y), r, r, r);
620       
621 #if 0
622       spec.Form("<area shape='rect' alt='%d' title='%d' href='qa_%09d.html' "
623                 "coords='%d,%d,%d,%d'>", r, r, r, 
624                 UInt_t(cw*(x-tx)), 0, UInt_t(cw*(x+tx)), ch);
625 #endif
626       areas->Add(area);
627     }
628   }
629   /** 
630    * Output list of runs 
631    * 
632    * @param o Output stream
633    */
634   void WriteRuns(std::ostream& o)
635   {
636     o << "<div class='jobid'><!--JOBID--></div>\n"
637       << "<div class='runs'>\n"
638       << "  Runs:<br> ";
639     for (Int_t i = 0; i < fRuns.GetSize(); i++) {
640       o << "<a href='" << Form("%09d", fRuns[i]) << "/index.html'>"
641         << fRuns[i] << "</a> " << std::flush;
642     }
643     o << "\n"
644       << "</div>" << std::endl;
645   }
646   /** 
647    * Write page footer 
648    * 
649    */
650   void WriteFooter()
651   {
652     WriteRuns(*fHtml);
653     QABase::WriteFooter();
654   }
655   /** 
656    * Write out image footer 
657    * 
658    * @param o Output stream
659    * @param pngName Name of the PNG file 
660    */
661   void WriteImageFooter(std::ostream& o, const char* pngName)
662   {
663     WriteRuns(o);
664     QABase::WriteImageFooter(o, pngName);
665   }
666   TGraph*       fNAccepted; // Graph of number of accepted events
667   TGraphErrors* fVz;        // Graph of mean vertex
668   UInt_t        fFirst;     // First run
669   UInt_t        fLast;      // Last run
670   TArrayI       fRuns;      // Seen runs 
671   TObjArray     fFiles;
672   Bool_t        fUseVar;    // Use variance rather than min/max 
673 };
674  
675 //
676 // EOF
677 //
678
679   
680