]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/qa/QAPlotter.C
5cbaa5262cadac83ad687074bce373a78315961b
[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   /** 
199    * Compatiblity constructor 
200    * 
201    * @param prodYear    Year
202    * @param prodLetter  Letter
203    * @param useVar      Use variance 
204    */
205   QAPlotter(Long_t prodYear, Char_t prodLetter, Bool_t useVar) 
206     : QABase("", (prodYear < 2000 ? 2000 : 0) + prodYear, 
207              Form("LHC%02d%c", int(prodYear % 100), prodLetter), "pass0"), 
208       fNAccepted(0),
209       fVz(0), 
210       fUseVar(useVar)
211   {
212     Info("QAPlotter", "Do we use variance? %s", fUseVar ? "yes" : "no");
213     fFMD1i = new Ring(1, 'I', useVar); 
214     fFMD2i = new Ring(2, 'I', useVar); 
215     fFMD2o = new Ring(2, 'O', useVar); 
216     fFMD3i = new Ring(3, 'I', useVar); 
217     fFMD3o = new Ring(3, 'O', useVar); 
218     fNAccepted = new TGraph;
219     fNAccepted->SetName("nAccepted");
220     fNAccepted->SetMarkerStyle(20);
221     fNAccepted->SetLineWidth(2);
222
223     fVz = new TGraphErrors;
224     fVz->SetName("vz");
225     fVz->SetMarkerStyle(20);
226     fVz->SetLineWidth(2);
227   }
228   /** 
229    * Constructor 
230    */
231   QAPlotter(const TString& dataType, 
232             Int_t          year, 
233             const TString& period, 
234             const TString& pass, 
235             Bool_t         useVar=true) 
236     : QABase(dataType, year, period, pass),
237       fNAccepted(0),
238       fVz(0),
239       fUseVar(useVar)
240   {
241     Info("QAPlotter", "Do we use variance? %s", fUseVar ? "yes" : "no");
242     fFMD1i = new Ring(1, 'I', useVar); 
243     fFMD2i = new Ring(2, 'I', useVar); 
244     fFMD2o = new Ring(2, 'O', useVar); 
245     fFMD3i = new Ring(3, 'I', useVar); 
246     fFMD3o = new Ring(3, 'O', useVar); 
247     fNAccepted = new TGraph;
248     fNAccepted->SetName("nAccepted");
249     fNAccepted->SetMarkerStyle(20);
250     fNAccepted->SetLineWidth(2);
251
252     fVz = new TGraphErrors;
253     fVz->SetName("vz");
254     fVz->SetMarkerStyle(20);
255     fVz->SetLineWidth(2);
256   }
257   /** 
258    * Add a file to be processed
259    * 
260    * @param filename Name of file 
261    */
262   void AddFile(const char* filename)
263   {
264     fFiles.Add(new TObjString(filename));
265   }
266   const char* GetUserInfo(TList* l, const char* name, const char* def) const
267   {
268     if (!l) {
269       Warning("GetUserInfo", "No user information list");
270       return def;
271     }
272     
273     TObject* o = l->FindObject(name);
274     if (!o) {
275       Warning("GetUserInfo", "User information %s not found", name);
276       l->ls();
277       return def;
278     }
279
280     Info("GetUserInfo", "Got user information %s=%s", name, o->GetTitle());
281     return o->GetTitle();
282   }
283   /** 
284    * Make a tree 
285    * 
286    * @param read If true, read from file 
287    * 
288    * @return True on success 
289    */
290   Bool_t MakeTree(bool read)
291   {
292     if (fFiles.GetEntriesFast() <= 0) return QABase::MakeTree(read);
293     if (fFiles.GetEntriesFast() == 1 && read) {
294       TFile* file = TFile::Open(fFiles.At(0)->GetName(), "READ");
295       if (file) {
296         fTree       = static_cast<TTree*>(file->Get("T"));
297         return true;
298       }
299     }
300     TChain* chain = new TChain("T", "T");
301     if (!chain->AddFileInfoList(&fFiles)) return false;
302     
303     fTree   = chain;
304     
305     return true;
306   }
307
308   /** 
309    * Run the job
310    * 
311    */
312   void Run()
313   {
314     Init(true);
315     if (!fTree) {
316       Error("Run", "No input tree");
317       return;
318     }
319     fFirst = 0xFFFFFFFF;
320     fLast  = 0;
321
322     UInt_t nEntries = fTree->GetEntries();
323     UInt_t j = 0;
324     fRuns.Set(nEntries);
325     Info("Run", "Got %d runs", nEntries);
326     TList* l = fTree->GetUserInfo();
327     fPeriod   = GetUserInfo(l, "period", "?");
328     fPass     = GetUserInfo(l, "pass",   "?");
329     fDataType = GetUserInfo(l, "type",   "?");
330       
331
332     for (UInt_t i = 0; i < nEntries; i++) {
333       fTree->GetEntry(i);
334
335       UInt_t run = fGlobal->runNo;
336       UInt_t nev = fGlobal->nAccepted;
337
338       fFirst = TMath::Min(run, fFirst);
339       fLast  = TMath::Max(run, fLast);
340       fRuns[i] = run;
341
342       Info("Run", "Got run %d with %d accepted events", run, nev);
343       fNAccepted->SetPoint(i, run, nev);
344       fVz->SetPoint(i, run, fGlobal->meanVz);
345       fVz->SetPointError(i, 0, fGlobal->sigmaVz);
346
347       if (nev <= 100) continue;
348       static_cast<Ring*>(fFMD1i)->Update(j, run);
349       static_cast<Ring*>(fFMD2i)->Update(j, run);
350       static_cast<Ring*>(fFMD2o)->Update(j, run);
351       static_cast<Ring*>(fFMD3i)->Update(j, run);
352       static_cast<Ring*>(fFMD3o)->Update(j, run);
353       j++;
354     }
355     
356     Plot();
357   }
358   /** 
359    * Plot results
360    * 
361    */
362   void Plot()
363   {
364     // fTeXName = Form("trend_%09d_%09d", fFirst, fLast);
365     TString title;
366     if (!fPeriod.IsNull() && !fPass.IsNull()) 
367       title.Form("QA trends for %s/%s runs %d --- %d", 
368                  fPeriod.Data(), fPass.Data(), fFirst, fLast);
369     else
370       title.Form("QA trends for runs %d --- %d", fFirst, fLast);
371     MakeCanvas(title);
372
373     CanvasTitle("# of accepted events");
374     fNAccepted->Draw("apl");
375     PutCanvasTitle("# of accepted events");
376     AddRuns(fNAccepted->GetHistogram(), "# of accepted events");
377
378     TLine* l = new TLine(fFirst, 100, fLast, 100);
379     l->SetLineColor(kRed+2);
380     l->SetLineStyle(2);
381     l->Draw();
382     PrintCanvas("nAccepted");
383
384     CanvasTitle("#LTv_{z}#GT");
385     fVz->Draw("apl");
386     PutCanvasTitle("Mean z coordinate of interaction point");
387     AddRuns(fVz->GetHistogram(), "#LTv_{z}#GT");
388     PrintCanvas("vz");
389
390     TMultiGraph* chi2           = new TMultiGraph;              
391     TMultiGraph* c              = new TMultiGraph;              
392     TMultiGraph* delta          = new TMultiGraph;              
393     TMultiGraph* xi             = new TMultiGraph;              
394     TMultiGraph* sigma          = new TMultiGraph;              
395     TMultiGraph* low            = new TMultiGraph;              
396     TMultiGraph* singles        = new TMultiGraph;      
397     TMultiGraph* loss           = new TMultiGraph;
398     TMultiGraph* occ    = new TMultiGraph;      
399     TMultiGraph* beta           = new TMultiGraph;              
400     chi2        ->SetName("chi2");             
401     c           ->SetName("c");                
402     delta       ->SetName("delta");            
403     xi          ->SetName("xi");               
404     sigma       ->SetName("sigma");            
405     low         ->SetName("low");              
406     singles     ->SetName("singles");    
407     loss        ->SetName("loss");
408     beta        ->SetName("beta");          
409     occ ->SetName("occupancy");
410
411
412     AddToMulti(fFMD1i,chi2, c, delta, xi, sigma, low, singles, loss, beta, occ);
413     AddToMulti(fFMD2i,chi2, c, delta, xi, sigma, low, singles, loss, beta, occ);
414     AddToMulti(fFMD2o,chi2, c, delta, xi, sigma, low, singles, loss, beta, occ);
415     AddToMulti(fFMD3i,chi2, c, delta, xi, sigma, low, singles, loss, beta, occ);
416     AddToMulti(fFMD3o,chi2, c, delta, xi, sigma, low, singles, loss, beta, occ);
417
418     PlotMulti(chi2,     "#LT#chi^{2}/#nu#GT from #Delta fits", true);
419     PlotMulti(c,        "#LTc#GT from #Delta fits");
420     PlotMulti(delta,    "#LT#Delta_{p}#GT from #Delta fits");
421     PlotMulti(xi,       "#LT#xi#GT from #Delta fits");
422     PlotMulti(sigma,    "#LT#sigma#GT from #Delta fits");
423     PlotMulti(low,      "Bins with too low statistics");
424     PlotMulti(singles,  "Fraction of single hits");
425     PlotMulti(loss,     "% of hits 'lost' due to merging+cuts");
426     PlotMulti(occ,      "#LTOccupancy#GT [%]", true);
427     PlotMulti(beta,     "Correlation of methods");
428
429     fStore->cd();
430     fNAccepted->Write();
431     chi2->Write();
432     c->Write();
433     delta->Write();
434     sigma->Write();
435     low->Write();
436     singles->Write();
437     loss->Write();
438     occ->Write();
439     beta->Write();
440
441     std::ofstream doc(".doc");
442     doc << fPeriod << " " << fPass << " ("
443         << fDataType << ")" << std::endl;
444     doc.close();
445
446     Close(false); // Do not delete PNGs
447   }
448   /** 
449    * Add graphs from a ring to the multi graphs
450    * 
451    * @param qr        Ring
452    * @param chi2      ELoss reduced chi-square     
453    * @param c         ELoss constant               
454    * @param delta     ELoss MPV            
455    * @param xi        ELoss Landau width           
456    * @param sigma     ELoss Gaus width             
457    * @param low       bins with low statistics     
458    * @param singles   fraction of singles          
459    * @param loss      'lost' data                  
460    * @param beta      Poisson vs ELoss correlation
461    * @param occupancy mean occupancy              
462    */
463   void AddToMulti(QARing* qr, 
464                   TMultiGraph* chi2,
465                   TMultiGraph* c,
466                   TMultiGraph* delta,
467                   TMultiGraph* xi,
468                   TMultiGraph* sigma,
469                   TMultiGraph* low,
470                   TMultiGraph* singles,
471                   TMultiGraph* loss,
472                   TMultiGraph* beta,
473                   TMultiGraph* occupancy)
474   {
475     Ring* r = static_cast<Ring*>(qr);
476     chi2        ->Add(r->fGChi2);
477     c           ->Add(r->fGC);
478     delta       ->Add(r->fGDelta);
479     xi          ->Add(r->fGXi);
480     sigma       ->Add(r->fGSigma);
481     low         ->Add(r->fGLow);
482     singles     ->Add(r->fGSingles);
483     loss        ->Add(r->fGLoss);
484     occupancy   ->Add(r->fGOccupancy);
485     beta        ->Add(r->fGBeta);
486   }
487   /** 
488    * Plot a multi-graph
489    * 
490    * @param mg     Multi graph
491    * @param title  Title
492    * @param logy   If true, make @f$\log@f$ scale 
493    */
494   void PlotMulti(TMultiGraph* mg, const char* title, Bool_t logy=false)
495   {
496     CanvasTitle(title);
497     // fCanvas->SetBottomMargin(.15);
498     fCanvas->SetLeftMargin(.08);
499     fCanvas->SetTopMargin(.1);
500     fCanvas->SetLogy(logy);
501     // fCanvas->SetRightMargin(.2);
502     mg->Draw("apl");
503     
504     TH1*     h   = mg->GetHistogram();
505     Double_t max = h->GetMaximum();
506     Double_t min = h->GetMinimum();
507     if (h->GetMinimum() == 0) {
508       min = min - .1*(max-min);
509       h->SetMinimum(min);
510     }
511     Int_t    x1  = h->GetXaxis()->GetXmin();
512     Int_t    x2  = h->GetXaxis()->GetXmax();
513     
514     TLegend* l = new TLegend(.1, .91, .97, .95);
515     l->SetNColumns(5);
516     l->SetFillColor(0);
517     l->SetFillStyle(0);
518     l->SetBorderSize(0);
519     l->SetTextFont(42);
520     TIter next(mg->GetListOfGraphs());
521     mg->GetListOfGraphs();
522     TGraph* g = 0;
523
524     // Get the runs we have here 
525     TArrayI runs(fRuns.GetSize());
526     runs.Reset(INT_MAX);
527     Int_t   j = 0;
528     while ((g = static_cast<TGraph*>(next()))) { 
529       l->AddEntry(g, g->GetTitle(), "lp");
530       Double_t* xs = g->GetX();
531       Int_t     n  = g->GetN();
532       for (Int_t i = 0; i < n; i++) {
533         if (FindRun(runs, Int_t(xs[i])) >= 0) continue;
534         runs.SetAt(xs[i], j++);
535       }
536       Double_t ymean = g->GetMean(2);
537       Double_t xh    = x2 - .03 * Double_t(x2-x1); 
538       TLine* lm = new TLine(x1, ymean, xh, ymean);
539       lm->SetLineColor(g->GetLineColor());
540       lm->SetLineStyle(2);
541       lm->SetLineWidth(1);
542       lm->Draw();
543       TLatex* al = new TLatex(xh, ymean, g->GetTitle());
544       al->SetTextColor(g->GetLineColor());
545       al->SetTextFont(42);
546       al->SetTextSize(.02);
547       al->SetTextAlign(12);
548       al->Draw();
549     }
550     l->Draw();
551
552     TList areas;
553     areas.SetOwner();
554     AddRuns(h, title, &runs, &areas);
555
556     PrintCanvas(mg->GetName(), &areas);
557   }
558   /** 
559    * Find a run 
560    * 
561    * @param runs List of runs 
562    * @param run  Run to find
563    * 
564    * @return Index of run in run list
565    */
566   Int_t FindRun(const TArrayI& runs, Int_t run) 
567   {
568     std::sort(&(runs.fArray[0]), &(runs.fArray[runs.GetSize()]));
569     Int_t idx = TMath::BinarySearch(runs.GetSize(), runs.fArray, run);
570     if (idx >= runs.GetSize() || idx < 0 || runs[idx] != run) return -1;
571     return idx;
572   }
573   /** 
574    * Add run labels at appropriate places on the plot
575    * 
576    * @param h      Frame histogram
577    * @param title  Title 
578    * @param runs   List of runs, if any
579    * @param areas  Other areas 
580    */
581   void AddRuns(TH1* h, const char* title, TArrayI* runs=0, 
582                TList* areas=0)
583   {
584     h->GetXaxis()->SetNoExponent();
585     // h->GetXaxis()->SetTitleOffset(1);
586     TString ytitle(title);
587     if (fUseVar) ytitle.Append(" (errors: variance)");
588     else         ytitle.Append(" (errors: min/max)");
589     h->SetYTitle(ytitle.Data());
590     h->SetXTitle("Run #");
591     // Info("AddRuns", "%s: %s vs %s", h->GetName(), 
592     //      h->GetXaxis()->GetTitle(), h->GetYaxis()->GetTitle());
593   
594     Int_t    r1  = h->GetXaxis()->GetXmin();
595     Int_t    r2  = h->GetXaxis()->GetXmax();
596     Double_t lx  = 0;
597     Double_t tx  = .025; // (r2 - r1) / 18;
598     Double_t wx  = 1 - fCanvas->GetLeftMargin() - fCanvas->GetRightMargin();
599     Double_t dy  = .025;
600     Double_t y   = fCanvas->GetBottomMargin()+dy;
601     for (Int_t i = 0; i < fRuns.GetSize(); i++) {
602       Int_t    r = fRuns[i];
603       Double_t x = fCanvas->GetLeftMargin() + wx*Double_t(r-r1)/(r2-r1);
604
605       // Skip runs out of range 
606       if (r < r1 || r > r2) continue;
607
608       // Skip runs not in the graphs 
609       if (runs) {
610         if (FindRun(*runs, r) < 0) { 
611           continue;
612         }
613       }
614           
615
616       if (TMath::Abs(x - lx) < tx) y += dy;
617       else                         y =  fCanvas->GetBottomMargin() + dy;
618
619
620       // Info(h->GetName(), "%6d (x,y)=(%12f,%12f) |lx-x|=|%f-x|=%f", 
621       //      r, x, y, lx, TMath::Abs(lx-x));
622       lx = x;
623       
624       Double_t* xa = fNAccepted->GetX();
625       Int_t     na = fNAccepted->GetN();
626       Int_t idx = TMath::BinarySearch(na, xa, Double_t(r));
627       Color_t color = kBlue+3;
628       if (idx >= 0  && idx < na  && r == xa[idx] && 
629           fNAccepted->GetY()[idx] < 10000) 
630         color = kRed+3;
631
632       TLatex* ll = new TLatex(x, y, Form("%d", r));
633       ll->SetNDC();
634       ll->SetTextAlign(21);
635       ll->SetTextSize(0.02);
636       ll->SetTextColor(color);
637       ll->SetTextFont(42);
638       // ll->SetTextAngle(90);
639       ll->Draw();
640       TLine* tl = new TLine(x, y, x, 1-fCanvas->GetTopMargin());
641       tl->SetBit(TLine::kLineNDC);
642       tl->SetLineStyle(3);
643       tl->SetLineColor(color);
644       tl->Draw();
645
646       if (!areas) continue;
647       
648       TObjString* area = new TObjString;
649       TString&    spec = area->String();
650       spec.Form("<span style=\"left: %d%%; bottom: %d%%;\" "
651                 "onClick='window.location=\"%09d/index.html\"' "
652                 "onMouseOver='this.style.cursor=\"pointer\"' "
653                 "alt=\"%d\""
654                 ">%d</span>",
655                 UInt_t(100*x)-2, UInt_t(100*y), r, r, r);
656       
657 #if 0
658       spec.Form("<area shape='rect' alt='%d' title='%d' href='qa_%09d.html' "
659                 "coords='%d,%d,%d,%d'>", r, r, r, 
660                 UInt_t(cw*(x-tx)), 0, UInt_t(cw*(x+tx)), ch);
661 #endif
662       areas->Add(area);
663     }
664   }
665   /** 
666    * Output list of runs 
667    * 
668    * @param o Output stream
669    */
670   void WriteRuns(std::ostream& o)
671   {
672     o << "<div class='jobid'><!--JOBID--></div>\n"
673       << "<div class='runs'>\n"
674       << "  Runs:<br> ";
675     for (Int_t i = 0; i < fRuns.GetSize(); i++) {
676       o << "<a href='" << Form("%09d", fRuns[i]) << "/index.html'>"
677         << fRuns[i] << "</a> " << std::flush;
678     }
679     o << "\n"
680       << "</div>" << std::endl;
681   }
682   /** 
683    * Write page footer 
684    * 
685    */
686   void WriteFooter()
687   {
688     WriteRuns(*fHtml);
689     QABase::WriteFooter();
690   }
691   /** 
692    * Write out image footer 
693    * 
694    * @param o Output stream
695    * @param pngName Name of the PNG file 
696    */
697   void WriteImageFooter(std::ostream& o, const char* pngName)
698   {
699     WriteRuns(o);
700     QABase::WriteImageFooter(o, pngName);
701   }
702   TGraph*       fNAccepted; // Graph of number of accepted events
703   TGraphErrors* fVz;        // Graph of mean vertex
704   UInt_t        fFirst;     // First run
705   UInt_t        fLast;      // Last run
706   TArrayI       fRuns;      // Seen runs 
707   TObjArray     fFiles;
708   Bool_t        fUseVar;    // Use variance rather than min/max 
709 };
710  
711 //
712 // EOF
713 //
714
715   
716