ROOT logo
//____________________________________________________________________
//
// $Id$
//
// Script that contains a class to draw hits, using the
// AliFMDInputHits class in the util library. 
//
// It draws the energy loss versus the p/(mq^2).  It can be overlayed
// with the Bethe-Bloc curve to show how the simulation behaves
// relative to the expected. 
//
// Use the script `Compile.C' to compile this class using ACLic. 
//
#include <TH2D.h>
#include <AliFMDHit.h>
#include <AliFMDInput.h>
#include <AliStack.h>
#include <iostream>
#include <TStyle.h>
#include <TArrayF.h>
#include <TParticle.h>
#include <TCanvas.h>
#include <TGraphErrors.h>
#include <TDatabasePDG.h>
#include <TParticlePDG.h>
#include <TLegend.h>
#include <TArrow.h>
#include <TLatex.h>
#include <TF1.h>

/** @class DrawHits
    @brief Draw hit energy loss
    @code 
    Root> .L Compile.C
    Root> Compile("DrawHits.C")
    Root> DrawHits c
    Root> c.Run();
    @endcode
    @ingroup FMD_script
 */
class DrawHits : public AliFMDInput
{
private:
  TH2D* fElossVsPMQ; // Histogram 
  TH1D* fEloss;
  TH1D* fBetaGamma;
  TParticlePDG* fPdg;
  const Double_t fRho;
  const Double_t fBetaGammaMip;
public:
  //__________________________________________________________________
  DrawHits(const char* pdgName="pi+",
	   Int_t m=500, Double_t emin=1, Double_t emax=1000, 
	   Int_t n=900, Double_t tmin=1e-2, Double_t tmax=1e3) 
    : AliFMDInput("galice.root"),
      fElossVsPMQ(0),
      fEloss(0),
      fBetaGamma(0),
      fPdg(0), 
      fRho(2.33), // 2.33),
      fBetaGammaMip(3.4601)
  { 
    AddLoad(kKinematics);
    AddLoad(kHits);
    TDatabasePDG* pdgDB = TDatabasePDG::Instance();
    fPdg                = pdgDB->GetParticle(pdgName);
    if (!fPdg) Warning("DrawHits", "Particle %s not found", pdgName);

    TArrayF tkine(MakeLogScale(n,    tmin, tmax));
    TArrayF betag(MakeLogScale(n/10, tmin, tmax));
    TArrayF eloss(MakeLogScale(m,    emin, emax));
    TString name("elossVsPMQ");
    TString title(Form("#Delta E/#Delta x / q^{2} vs. p/m, %s", 
		       (pdgName ? pdgName : "")));
    fElossVsPMQ = new TH2D(name.Data(), title.Data(), 
			   tkine.fN-1, tkine.fArray, 
			   eloss.fN-1, eloss.fArray);
    fElossVsPMQ->SetXTitle("p/(mq^{2})=#beta#gamma/q^{2}");
    fElossVsPMQ->SetYTitle("#Delta E/#Delta x / q^{2} [MeV/cm]");
    fElossVsPMQ->Sumw2();
    fEloss = new TH1D("eloss", "#Delta E/#Delta x / q^{2}", 
		      eloss.fN-1, eloss.fArray);
    fEloss->SetFillColor(kRed);
    fEloss->SetFillStyle(3001);
    fEloss->SetXTitle("#Delta E/#Delta x / q^{2} [MeV/cm]");
    fEloss->Sumw2();

    fBetaGamma = new TH1D("betaGamma", "#beta#gamma of particles", 
			  betag.fN-1, betag.fArray);
    fBetaGamma->SetXTitle("#beta#gamma");
    fBetaGamma->SetFillColor(kBlue+1);
    fBetaGamma->SetFillStyle(3001);
  }
  //__________________________________________________________________
  Bool_t ProcessHit(AliFMDHit* hit, TParticle* p) 
  {
    if (!hit) {
      std::cout << "No hit" << std::endl;
      return kFALSE;
    }

    if (!p) {
      std::cout << "No track" << std::endl;
      return kFALSE;
    }
    // if (!p->IsPrimary()) return kTRUE;
    if (hit->IsStop()) return kTRUE;
    if (hit->Length() == 0) { 
      Warning("ProcessHit", "Hit in %s has 0 length", hit->GetName());
      return kTRUE;
    }
    
    Float_t q = hit->Q() / 3.;
    Float_t m = hit->M();
    if (m == 0 || q == 0) return kTRUE;

    TLorentzVector pp;
    p->Momentum(pp);
    Double_t betagamma = 0;
    Info("ProcessHit", "%s (%s) beta=%f", p->GetPDG()->GetName(), 
	 fStack->IsPhysicalPrimary(hit->Track()) ? "primary" : "secondary", 
	 pp.Beta());
    if (pp.Beta() <= 1 && pp.Beta() >= 0) 
      betagamma = pp.Beta() * pp.Gamma();
    fBetaGamma->Fill(betagamma);
#if 0
    if (betagamma < 10) { 
      Info("ProcessHit", "%s (%s) beta=%f gamma=%f beta*gamma=%f", 
	   p->GetPDG()->GetName(), 
	   fStack->IsPhysicalPrimary(hit->Track()) ? "primary" : 
	   "secondary", 
	   pp.Beta(), pp.Gamma(), betagamma);
      return kTRUE;
    }
#endif

    Float_t x = hit->P();
    Float_t y = hit->Edep()/hit->Length();

    x /= hit->M();
    // y /= q * q;
    fElossVsPMQ->Fill(x, y);
    fEloss->Fill(y);
    // fNHits++;
    return kTRUE;
  }
  //__________________________________________________________________
  void ShowFit(Double_t x1, Double_t y1, const char* title, 
	       TF1* f, Double_t dx=0, Double_t dy=0.05)
  {
    Double_t x = x1, y = y1;
    TLatex* latex = new TLatex(x, y, title);
    latex->SetTextFont(132);
    latex->SetTextSize(0.8*dy);
    latex->SetNDC();
    latex->Draw();
    x -= dx;
    y -= dy;
    const Double_t eqDx=0.1;
    Double_t chi2 = f->GetChisquare();
    Int_t    ndf  = f->GetNDF();
    Double_t prob = f->GetProb();
    latex->DrawLatex(x, y, "#chi^{2}/NDF");
    latex->DrawLatex(x+eqDx, y, Form("= %7.4f/%3d=%5.2f (%3d%%)", 
				     chi2, ndf, chi2/ndf, int(100*prob)));
    Int_t     n = f->GetNpar();
    Double_t* p = f->GetParameters();
    Double_t* e = f->GetParErrors();
    for (int i = 0; i < n; i++) { 
      x -= dx;
      y -= dy;
      latex->DrawLatex(x, y, f->GetParName(i));
      latex->DrawLatex(x+eqDx, y, Form("= %7.4f", p[i]));
      latex->DrawLatex(x+2*eqDx, y, Form("#pm %7.4f", e[i]));
    }
  }      
  //__________________________________________________________________
  Bool_t Finish()
  {
    gStyle->SetPalette(1);
    // gStyle->SetOptTitle(0);
    gStyle->SetTitleBorderSize(1);
    gStyle->SetTitleFillColor(0);
    gStyle->SetCanvasColor(0);
    gStyle->SetCanvasColor(0);
    gStyle->SetCanvasBorderSize(0);
    gStyle->SetPadColor(0);
    gStyle->SetPadBorderSize(0);
    gStyle->SetOptStat(0);
    TCanvas* c = new TCanvas("elossVsP", "Energy loss versus momentum", 
			     1200, 800);
    c->SetLogy();
    c->SetLogx();

    TString title(Form("%s, %d events", fElossVsPMQ->GetTitle(), fEventCount));
    fElossVsPMQ->SetTitle(title.Data());
    fElossVsPMQ->SetStats(kFALSE);
    fElossVsPMQ->Draw("AXIS");
    fElossVsPMQ->Draw("ACOLZ same");
    TGraph*  mate    = FromGFMATE();
    TGraph*  bb      = FromRPPFull();
    TGraph*  nodelta = FromRPPNoDelta();
    TGraph*  norad   = FromRPPNoRad();
    TGraph*  mean    = FromRPPMean();
    // fElossVsPMQ->Draw("ACOLZ same");
    TLegend* l       = new TLegend(.5, .6, .89, .89);
    // l->AddEntry(fElossVsPMQ, fElossVsPMQ->GetTitle(), "pf");
    l->SetFillColor(0);
    l->AddEntry(mate,        mate->GetTitle(),    "lf");
    l->AddEntry(bb,          bb->GetTitle(),      "l");
    l->AddEntry(nodelta,     nodelta->GetTitle(), "l");
    l->AddEntry(norad,       norad->GetTitle(),   "l");
    l->AddEntry(mean,        mean->GetTitle(),    "l");
    l->Draw("same");
    Double_t min = fElossVsPMQ->GetYaxis()->GetFirst();
    TArrow* a = new TArrow(fBetaGammaMip, min, fBetaGammaMip, 35*min, 03, "<|");
    a->SetAngle(30);
    a->Draw("same");
    TLatex* t = new TLatex(fBetaGammaMip, 40*min, "Minimum Ionising");
    t->SetTextSize(0.04);
    t->SetTextAlign(21);
    t->Draw("same");
    c->Modified();
    c->Update();
    c->cd();
    c->SaveAs("eloss_bethe.png");

    c = new TCanvas("cEloss", "Energy loss per unit material",
		    1200, 800);
    c->SetRightMargin(0.05);
    c->SetTopMargin(0.05);
    c->SetLeftMargin(0.05);
    fEloss->SetStats(kFALSE);
    // c->SetLogx();
    TF1* land     = new TF1("land", "landau");
    land->SetParNames("A", "MPV", "width");
    land->SetLineWidth(2);
    land->SetLineColor(kGreen+1);

    TF1* landgaus = new TF1("landgaus", "gaus(0)+landau(3)");
    landgaus->SetParNames("A", "#mu", "#sigma", "B", "MPV", "width");
    landgaus->SetLineWidth(2);
    landgaus->SetLineColor(kMagenta+1);
    TGraph*  corr  = GetCorr();
    TGraph*  resp  = GetResp();
    if (fEloss->GetEntries() != 0) { 
      c->SetLogy();
      fEloss->Scale(1. / fEloss->GetEntries());
      fEloss->GetXaxis()->SetRangeUser(1, 10);
      fEloss->Fit(land, "+", "", 2, 10);
      landgaus->SetParameters(land->GetParameter(0) / 100, 
			      land->GetParameter(1), 
			      land->GetParameter(2), 
			      land->GetParameter(0),
			      land->GetParameter(1), 
			      land->GetParameter(2));
      fEloss->Fit(landgaus, "+", "", 1, 10);
      fEloss->DrawCopy("HIST SAME");
    }
    
    // fEloss->DrawCopy("E SAME");
    // land->Draw("same");
    // landgaus->Draw("same");
    Double_t max = fEloss->GetMaximum();
    Double_t* x   = resp->GetX();
    Double_t* y   = resp->GetY();
    TGraph*   g   = new TGraph(resp->GetN());
    g->SetName(Form("%sCorr", resp->GetName()));
    g->SetTitle(resp->GetTitle());
    g->SetLineStyle(resp->GetLineStyle());
    g->SetLineColor(resp->GetLineColor());
    g->SetLineWidth(resp->GetLineWidth());
    Double_t  xs2 = corr->Eval(fBetaGammaMip);
    Double_t xss   = 1.1;
    Double_t  xs  = fRho * xss;
    std::cout << "Correction factor: " << xs2 << std::endl;
    for (Int_t i = 0; i < g->GetN(); i++) 
      g->SetPoint(i, x[i] * xs, y[i] * max);
    g->Draw("C same");
    
    l = new TLegend(.05, .6, .4, .95);
    l->SetFillColor(0);
    l->SetBorderSize(1);
    l->AddEntry(fEloss, fEloss->GetTitle(), "lf");
    if (land) 
      l->AddEntry(land,   Form("Landau fit\t- #chi^{2}/NDF=%7.5f", 
			       land->GetChisquare()/land->GetNDF()), "l");
    if (landgaus) 
      l->AddEntry(landgaus,   
		  Form("Landau+Gauss fit\t- #chi^{2}/NDF=%7.5f", 
		       landgaus->GetChisquare()/landgaus->GetNDF()), "l");
    l->AddEntry(resp,   Form("f(%s#Delta/x) 320#mum Si [RPP fig 27.7]",
			     xss != 1 ? Form("%4.2f#times", xss) : ""), 
		"l");
    l->Draw("same");
    
    fEloss->GetYaxis()->SetRangeUser(1e-4, 100);
    ShowFit(0.45,.9, "Landau+Gaus", landgaus, 0, 0.02);
    ShowFit(0.70,.9, "Landau", land, 0, 0.02);

    c->Modified();
    c->Update();
    c->cd();
    c->SaveAs("eloss_landau.png");


    c = new TCanvas("cBetaGamma", "beta gamma", 1200, 800);
    c->SetLogx();
    fBetaGamma->Draw();
    Int_t mipbin = fBetaGamma->FindBin(fBetaGammaMip) + 1;
    Int_t maxbin = fBetaGamma->GetNbinsX();
    Int_t total  = fBetaGamma->Integral();
    Int_t over   = fBetaGamma->Integral(mipbin,maxbin);
    TH1*  res    = (TH1*)fBetaGamma->Clone("overMip");
    res->SetFillColor(kRed+1);
    for (int i = 0; i < mipbin; i++) 
      res->SetBinContent(i, 0);
    res->Draw("same");
    std::cout << "Percentage over MIP : " << float(over) / total << std::endl;

    Double_t yy = fBetaGamma->GetBinContent(mipbin) * 1.1; 
    a = new TArrow(fBetaGammaMip, 0, fBetaGammaMip, yy, 3, "<|");
    a->SetAngle(30);
    a->Draw("same");
    t = new TLatex(fBetaGammaMip, yy, "Minimum Ionising");
    t->SetTextSize(0.04);
    t->SetTextAlign(21);
    t->Draw("same");
    c->Modified();
    c->Update();
    c->cd();
    c->SaveAs("eloss_betagamma.png");
    
    return kTRUE;
  }

  /** Scale a graph by density (multiply) and mass (divide). 
      @param graph Graph to scale 
      @param density If @c true, scale by the Si density
      ($\rho=2.33$/cm^3$).  The y axis is assumed to have units of
      $MeVg^{-1}cm^2$. 
      @param mass Mass to scale with. The x axis is assumed to be the
      kinetic energy of the particles in units of $GeV$.  */
  void ScaleGraph(TGraph* graph, bool density=true, double mass=1) 
  {
      Double_t*      x   = graph->GetX();
      Double_t*      y   = graph->GetY();
      const Double_t rho = (density ? fRho : 1);
      for (Int_t i = 0; i < graph->GetN(); i++) 
	graph->SetPoint(i, x[i] / mass, y[i] * rho); 
  }    
  /** Draw pure Bethe-Bloc from Review of Particle Physics, fig. 27.1 
      @return TGraph object */ 
  TGraph* FromRPPFull() 
  {
    static TGraph* graph = 0;
    if (!graph) { 
      graph = new TGraph(20);
      graph->GetHistogram()->SetXTitle("#beta#gamma");
      graph->GetHistogram()->SetYTitle("#Delta E/#Delta x [MeV/cm]");
      graph->SetFillColor(0);
      graph->SetLineColor(kRed+1);
      graph->SetLineStyle(2);
      graph->SetLineWidth(2);
      graph->SetName("full_stop");
      graph->SetTitle("Stopping (MeVcm^{2}/g) [RPP fig 27.1]");
      graph->SetPoint(0,0.001461622,40.17542);
      graph->SetPoint(1,0.003775053,91.28429);
      graph->SetPoint(2,0.01178769,202.7359);
      graph->SetPoint(3,0.01722915,212.1938);
      graph->SetPoint(4,0.03162278,172.8318);
      graph->SetPoint(5,0.06028646,91.28429);
      graph->SetPoint(6,0.09506529,51.62633);
      graph->SetPoint(7,0.433873,5.281682);
      graph->SetPoint(8,1.255744,1.808947);
      graph->SetPoint(9,2.393982,1.440177);
      graph->SetPoint(10,3.499097,1.407715);
      graph->SetPoint(11,10.92601,1.542122);
      graph->SetPoint(12,60.28646,1.85066);
      graph->SetPoint(13,236.3885,2.121938);
      graph->SetPoint(14,468.0903,2.324538);
      graph->SetPoint(15,1208.976,2.987085);
      graph->SetPoint(16,6670.768,7.961412);
      graph->SetPoint(17,23341.67,24.3298);
      graph->SetPoint(18,110651.2,104.6651);
      graph->SetPoint(19,264896.9,260.5203);
      ScaleGraph(graph);
    }
    graph->Draw("C same");
    return graph;
  }

  /** Draw pure Bethe-Bloc from Review of Particle Physics, fig. 27.1,
      but without delta electrons  
      @return TGraph object */ 
  TGraph* FromRPPNoDelta() 
  {
    static TGraph* graph = 0;
    if (!graph) { 
      graph = new TGraph(20);
      graph->SetName("stop_nodelta");
      graph->SetTitle("Stopping w/o #delta's [RPP fig 27.1]");
      graph->GetHistogram()->SetYTitle("(MeVcm^{2}/g)");
      graph->GetHistogram()->SetXTitle("#beta#gamma");
      graph->SetFillColor(0);
      graph->SetLineColor(kGreen+1);
      graph->SetLineStyle(2);
      graph->SetLineWidth(2);
      graph->SetPoint(0,0.001461622,40.17542);
      graph->SetPoint(1,0.003775053,91.28429);
      graph->SetPoint(2,0.01178769,202.7359);
      graph->SetPoint(3,0.01722915,212.1938);
      graph->SetPoint(4,0.03162278,172.8318);
      graph->SetPoint(5,0.06028646,91.28429);
      graph->SetPoint(6,0.09506529,51.62633);
      graph->SetPoint(7,0.433873,5.281682);
      graph->SetPoint(8,1.255744,1.808947);
      graph->SetPoint(9,2.304822,1.473387);
      graph->SetPoint(10,3.921088,1.473387);
      graph->SetPoint(11,8.064796,1.614064);
      graph->SetPoint(12,26.15667,1.936996);
      graph->SetPoint(13,264.8969,2.489084);
      graph->SetPoint(14,544.8334,2.665278);
      graph->SetPoint(15,1163.949,2.853945);
      graph->SetPoint(16,5312.204,3.19853);
      graph->SetPoint(17,15374.93,3.424944);
      graph->SetPoint(18,49865.73,3.667384);
      graph->SetPoint(19,634158.5,4.110185);
      ScaleGraph(graph);
    }
    graph->Draw("C same");
    return graph;
  }

  /** Draw pure Bethe-Bloc from Review of Particle Physics, fig. 27.1,
      but without delta electrons  
      @return TGraph object */ 
  TGraph* FromRPPNoRad() 
  {
    static TGraph* graph = 0;
    if (!graph) { 
      graph = new TGraph(18);
      graph->SetName("norad_stop");
      graph->SetTitle("Stopping w/o radiative loss [RPP fig. 27.1]");
      graph->GetHistogram()->SetYTitle("(MeVcm^{2}/g)");
      graph->GetHistogram()->SetXTitle("#beta#gamma");
      graph->SetFillColor(0);
      graph->SetLineColor(kBlue+1);
      graph->SetLineWidth(2);
      graph->SetLineStyle(2);
      graph->SetPoint(0,0.001,24.3298);
      graph->SetPoint(1,0.003117649,74.35105);
      graph->SetPoint(2,0.008675042,172.8318);
      graph->SetPoint(3,0.01782497,212.1938);
      graph->SetPoint(4,0.02704573,189.3336);
      graph->SetPoint(5,0.07481082,70.29816);
      graph->SetPoint(6,0.3300035,8.524974);
      graph->SetPoint(7,0.819559,2.489084);
      graph->SetPoint(8,1.447084,1.651284);
      graph->SetPoint(9,2.555097,1.440177);
      graph->SetPoint(10,4.026598,1.407715);
      graph->SetPoint(11,32.38084,1.728318);
      graph->SetPoint(12,97.19733,1.893336);
      graph->SetPoint(13,1732.539,2.170869);
      graph->SetPoint(14,11098.58,2.324538);
      graph->SetPoint(15,32075.46,2.378141);
      graph->SetPoint(16,221655.8,2.546482);
      graph->SetPoint(17,593830.6,2.605203);
      ScaleGraph(graph);
    }
    graph->Draw("C same");
    return graph;
  }

  /** Draw pure Bethe-Bloc from Review of Particle Physics, fig. 27.6 
      @return TGraph object */ 
  TGraph* FromRPPMean() 
  {
    static TGraph* graph = 0;
    if (!graph) { 
      graph = new TGraph(12);
      graph->SetName("mean_eloss");
      graph->SetTitle("Mean #Delta E/#Delta x - "
		      "electronic only  [RPP fig. 27.6]");
      graph->GetHistogram()->SetYTitle("(MeVcm^{2}/g)");
      graph->GetHistogram()->SetXTitle("#mu E_{kin} (GeV)");
      graph->SetFillColor(0);
      graph->SetLineStyle(2);
      graph->SetLineWidth(2);
      graph->SetLineColor(kMagenta+1);
      graph->SetMarkerStyle(21);
      graph->SetMarkerSize(0.6);
      graph->SetPoint(0,0.1,1.346561);
      graph->SetPoint(1,0.1435819,1.230159);
      graph->SetPoint(2,0.2061576,1.156085);
      graph->SetPoint(3,0.3698076,1.124339);
      graph->SetPoint(4,0.4620113,1.124339);
      graph->SetPoint(5,0.8521452,1.145503);
      graph->SetPoint(6,1.909707,1.177249);
      graph->SetPoint(7,4.048096,1.198413);
      graph->SetPoint(8,12.66832,1.219577);
      graph->SetPoint(9,48.17031,1.230159);
      graph->SetPoint(10,285.8863,1.230159);
      graph->SetPoint(11,894.6674,1.230159);
      const Double_t m   = 0.10566; // Muon 
      ScaleGraph(graph, true, m);
    }
    graph->Draw("C same");
    return graph;
  }

  /** Draw energy loss as obtained from GEANT 3.21 GFMATE. 
      @return TGraph object */
  TGraph* FromGFMATE() 
  {
    static TGraphErrors* gre = 0;
    if (!gre) {
      gre = new TGraphErrors(91);
      gre->SetName("ELOSS");
      gre->SetTitle("Energy loss 300#mu Si [GFMATE]");
      gre->GetHistogram()->SetXTitle("#beta#gamma");
      gre->GetHistogram()->SetYTitle("#Delta E/#Delta x [MeV/cm]");
      gre->SetFillColor(kGray);
      gre->SetFillStyle(3001); // 0 /* 3002 */);
      gre->SetLineColor(kGray+1);
      gre->SetLineStyle(1);
      gre->SetLineWidth(2);
      gre->SetPoint(0,7.16486e-05,1218.84);
      gre->SetPoint(1,9.25378e-05,1221.38);
      gre->SetPoint(2,0.000119517,1180.12);
      gre->SetPoint(3,0.000154362,1100.31);
      gre->SetPoint(4,0.000199367,996.621);
      gre->SetPoint(5,0.000257492,886.005);
      gre->SetPoint(6,0.000332563,780.483);
      gre->SetPoint(7,0.000429522,684.927);
      gre->SetPoint(8,0.000554749,599.407);
      gre->SetPoint(9,0.000716486,522.375);
      gre->SetPoint(10,0.000925378,452.497);
      gre->SetPoint(11,0.00119517,389.101);
      gre->SetPoint(12,0.00154362,331.974);
      gre->SetPoint(13,0.00199367,280.969);
      gre->SetPoint(14,0.00257492,235.689);
      gre->SetPoint(15,0.00332564,196.156);
      gre->SetPoint(16,0.00429522,162.402);
      gre->SetPoint(17,0.00554749,133.87);
      gre->SetPoint(18,0.00716486,109.959);
      gre->SetPoint(19,0.00925378,90.2035);
      gre->SetPoint(20,0.0119517,74.1317);
      gre->SetPoint(21,0.0154362,60.8988);
      gre->SetPoint(22,0.0199367,49.9915);
      gre->SetPoint(23,0.0257492,40.9812);
      gre->SetPoint(24,0.0332564,33.5739);
      gre->SetPoint(25,0.0429522,27.5127);
      gre->SetPoint(26,0.0554749,22.5744);
      gre->SetPoint(27,0.0716486,18.5674);
      gre->SetPoint(28,0.0925378,15.3292);
      gre->SetPoint(29,0.119517,12.7231);
      gre->SetPoint(30,0.154362,10.6352);
      gre->SetPoint(31,0.199367,8.97115);
      gre->SetPoint(32,0.257492,7.65358);
      gre->SetPoint(33,0.332564,6.61909);
      gre->SetPoint(34,0.429522,5.79837);
      gre->SetPoint(35,0.554749,5.148);
      gre->SetPoint(36,0.716486,4.65024);
      gre->SetPoint(37,0.925378,4.27671);
      gre->SetPoint(38,1.19517,3.99831);
      gre->SetPoint(39,1.54362,3.79877);
      gre->SetPoint(40,1.99367,3.6629);
      gre->SetPoint(41,2.57492,3.57594);
      gre->SetPoint(42,3.32564,3.52565);
      gre->SetPoint(43,4.29522,3.50206);
      gre->SetPoint(44,5.54749,3.49715);
      gre->SetPoint(45,7.16486,3.50467);
      gre->SetPoint(46,9.25378,3.51988);
      gre->SetPoint(47,11.9517,3.53932);
      gre->SetPoint(48,15.4362,3.56054);
      gre->SetPoint(49,19.9367,3.58189);
      gre->SetPoint(50,25.7492,3.60231);
      gre->SetPoint(51,33.2564,3.62113);
      gre->SetPoint(52,42.9522,3.638);
      gre->SetPoint(53,55.4749,3.65275);
      gre->SetPoint(54,71.6486,3.66537);
      gre->SetPoint(55,92.5378,3.67586);
      gre->SetPoint(56,119.517,3.68433);
      gre->SetPoint(57,154.362,3.69105);
      gre->SetPoint(58,199.367,3.6962);
      gre->SetPoint(59,257.492,3.69997);
      gre->SetPoint(60,332.564,3.70257);
      gre->SetPoint(61,429.522,3.70421);
      gre->SetPoint(62,554.749,3.70511);
      gre->SetPoint(63,716.486,3.7055);
      gre->SetPoint(64,925.378,3.70559);
      gre->SetPoint(65,1195.17,3.70558);
      gre->SetPoint(66,1543.62,3.70557);
      gre->SetPoint(67,1993.67,3.70555);
      gre->SetPoint(68,2574.92,3.70553);
      gre->SetPoint(69,3325.64,3.70552);
      gre->SetPoint(70,4295.22,3.7055);
      gre->SetPoint(71,5547.49,3.70548);
      gre->SetPoint(72,7164.86,3.70547);
      gre->SetPoint(73,9253.78,3.70545);
      gre->SetPoint(74,11951.7,3.70544);
      gre->SetPoint(75,15436.2,3.70544);
      gre->SetPoint(76,19936.7,3.70544);
      gre->SetPoint(77,25749.2,3.70544);
      gre->SetPoint(78,33256.4,3.70544);
      gre->SetPoint(79,42952.2,3.70544);
      gre->SetPoint(80,55474.9,3.70544);
      gre->SetPoint(81,71648.6,3.70544);
      gre->SetPoint(82,92537.8,3.70544);
      gre->SetPoint(83,119517,3.70544);
      gre->SetPoint(84,154362,3.70544);
      gre->SetPoint(85,199367,3.70544);
      gre->SetPoint(86,257492,3.70544);
      gre->SetPoint(87,332563,3.70544);
      gre->SetPoint(88,429522,3.70544);
      gre->SetPoint(89,554749,3.70544);
      gre->SetPoint(90,716486,3.70544);
      // Double_t* x = gre->GetX();
      Double_t* y = gre->GetY();
      for (Int_t i = 0; i < gre->GetN(); i++) 
	gre->SetPointError(i, 0, 2 * 0.1 * y[i]); // ! 1 sigma
    }
    gre->Draw("c3 same");
    return gre;
  }

  /** Get the response functin @f$ f(\Delta_p/x)@f$ from Review of
      Particle Physics (fig. 27.7).  It is scaled to the value at
      MPV. */ 
  TGraph* GetResp()
  {
    static TGraph*  graph = 0;
    if (!graph) {
      graph = new TGraph;
      graph->SetName("si_resp");
      graph->SetTitle("f(#Delta/x) scaled to the MPV value ");
      graph->GetHistogram()->SetXTitle("#Delta/x (MeVcm^{2}/g)");
      graph->GetHistogram()->SetYTitle("f(#Delta/x)");
      graph->SetLineColor(kBlue+1);
      graph->SetLineWidth(2);
      graph->SetFillColor(kBlue+1);
      graph->SetMarkerStyle(21);
      graph->SetMarkerSize(0.6);
#if 0
      // Figure 27.7 or Review of Particle physics - Straggeling function in 
      // silicon of 500 MeV Pions, normalised to unity at the most probable 
      // value.   
      graph->SetPoint(0,0.808094,0.00377358);
      graph->SetPoint(1,0.860313,0.0566038);
      graph->SetPoint(2,0.891645,0.116981);
      graph->SetPoint(3,0.912533,0.181132);
      graph->SetPoint(4,0.928198,0.260377);
      graph->SetPoint(5,0.938642,0.320755);
      graph->SetPoint(6,0.954308,0.377358);
      graph->SetPoint(7,0.964752,0.433962);
      graph->SetPoint(8,0.975196,0.490566);
      graph->SetPoint(9,0.98564,0.550943);
      graph->SetPoint(10,0.996084,0.611321);
      graph->SetPoint(11,1.00653,0.667925);
      graph->SetPoint(12,1.02219,0.732075);
      graph->SetPoint(13,1.03264,0.784906);
      graph->SetPoint(14,1.0483,0.845283);
      graph->SetPoint(15,1.06397,0.901887);
      graph->SetPoint(16,1.09008,0.958491);
      graph->SetPoint(17,1.10574,0.984906);
      graph->SetPoint(18,1.13708,1);
      graph->SetPoint(19,1.13708,1);
      graph->SetPoint(20,1.15796,0.988679);
      graph->SetPoint(21,1.17363,0.966038);
      graph->SetPoint(22,1.19974,0.916981);
      graph->SetPoint(23,1.2154,0.89434);
      graph->SetPoint(24,1.23629,0.837736);
      graph->SetPoint(25,1.2624,0.784906);
      graph->SetPoint(26,1.28329,0.724528);
      graph->SetPoint(27,1.3094,0.664151);
      graph->SetPoint(28,1.32507,0.611321);
      graph->SetPoint(29,1.3564,0.550943);
      graph->SetPoint(30,1.41384,0.445283);
      graph->SetPoint(31,1.44517,0.392453);
      graph->SetPoint(32,1.48695,0.335849);
      graph->SetPoint(33,1.52872,0.286792);
      graph->SetPoint(34,1.58094,0.237736);
      graph->SetPoint(35,1.63838,0.196226);
      graph->SetPoint(36,1.68016,0.169811);
      graph->SetPoint(37,1.75326,0.135849);
      graph->SetPoint(38,1.81593,0.113208);
      graph->SetPoint(39,1.89426,0.0981132);
      graph->SetPoint(40,1.96214,0.0830189);
      graph->SetPoint(41,2.0718,0.0641509);
      graph->SetPoint(42,2.19191,0.0490566);
      graph->SetPoint(43,2.31723,0.0415094);
      graph->SetPoint(44,2.453,0.0301887);
      graph->SetPoint(45,2.53133,0.0264151);
      graph->SetPoint(46,2.57833,0.0264151);
#else
      graph->SetPoint(0,0.8115124,0.009771987);
      graph->SetPoint(1,0.9198646,0.228013);
      graph->SetPoint(2,0.996614,0.5895765);
      graph->SetPoint(3,1.041761,0.8241042);
      graph->SetPoint(4,1.059819,0.8794788);
      graph->SetPoint(5,1.077878,0.9348534);
      graph->SetPoint(6,1.100451,0.980456);
      graph->SetPoint(7,1.141084,0.9967427);
      graph->SetPoint(8,1.204289,0.9153094);
      graph->SetPoint(9,1.276524,0.742671);
      graph->SetPoint(10,1.402935,0.465798);
      graph->SetPoint(11,1.515801,0.3029316);
      graph->SetPoint(12,1.73702,0.1465798);
      graph->SetPoint(13,1.985327,0.08143322);
      graph->SetPoint(14,2.301354,0.04234528);
      graph->SetPoint(15,2.56772,0.02931596);
#endif
    }
    return graph;
  }

  /** Get the correction to Bethe-Bloc from Review of Particle Physics
      (fig 27.8). 
  */
  TGraph* GetCorr() 
  {
    static TGraph* graph = 0;
    if (!graph) {
      graph = new TGraph(14);
      graph->SetName("graph");
      graph->SetTitle("(#Delta_{p}/x)/(dE/dx)|_{mip} for 320#mu Si");
      graph->GetHistogram()->SetXTitle("#beta#gamma = p/m");
      graph->SetFillColor(1);
      graph->SetLineColor(7);
      graph->SetMarkerStyle(21);
      graph->SetMarkerSize(0.6);
      graph->SetPoint(0,1.196058,0.9944915);
      graph->SetPoint(1,1.28502,0.9411017);
      graph->SetPoint(2,1.484334,0.8559322);
      graph->SetPoint(3,1.984617,0.7491525);
      graph->SetPoint(4,2.658367,0.6983051);
      graph->SetPoint(5,3.780227,0.6779661);
      graph->SetPoint(6,4.997358,0.6741525);
      graph->SetPoint(7,8.611026,0.684322);
      graph->SetPoint(8,15.28296,0.6995763);
      graph->SetPoint(9,41.54516,0.7186441);
      graph->SetPoint(10,98.91461,0.7288136);
      graph->SetPoint(11,203.2734,0.7326271);
      graph->SetPoint(12,505.6421,0.7338983);
      graph->SetPoint(13,896.973,0.7338983);
    }
    return graph;
  }
  
  ClassDef(DrawHits,0);
};

//____________________________________________________________________
//
// EOF
//
 DrawHits.C:1
 DrawHits.C:2
 DrawHits.C:3
 DrawHits.C:4
 DrawHits.C:5
 DrawHits.C:6
 DrawHits.C:7
 DrawHits.C:8
 DrawHits.C:9
 DrawHits.C:10
 DrawHits.C:11
 DrawHits.C:12
 DrawHits.C:13
 DrawHits.C:14
 DrawHits.C:15
 DrawHits.C:16
 DrawHits.C:17
 DrawHits.C:18
 DrawHits.C:19
 DrawHits.C:20
 DrawHits.C:21
 DrawHits.C:22
 DrawHits.C:23
 DrawHits.C:24
 DrawHits.C:25
 DrawHits.C:26
 DrawHits.C:27
 DrawHits.C:28
 DrawHits.C:29
 DrawHits.C:30
 DrawHits.C:31
 DrawHits.C:32
 DrawHits.C:33
 DrawHits.C:34
 DrawHits.C:35
 DrawHits.C:36
 DrawHits.C:37
 DrawHits.C:38
 DrawHits.C:39
 DrawHits.C:40
 DrawHits.C:41
 DrawHits.C:42
 DrawHits.C:43
 DrawHits.C:44
 DrawHits.C:45
 DrawHits.C:46
 DrawHits.C:47
 DrawHits.C:48
 DrawHits.C:49
 DrawHits.C:50
 DrawHits.C:51
 DrawHits.C:52
 DrawHits.C:53
 DrawHits.C:54
 DrawHits.C:55
 DrawHits.C:56
 DrawHits.C:57
 DrawHits.C:58
 DrawHits.C:59
 DrawHits.C:60
 DrawHits.C:61
 DrawHits.C:62
 DrawHits.C:63
 DrawHits.C:64
 DrawHits.C:65
 DrawHits.C:66
 DrawHits.C:67
 DrawHits.C:68
 DrawHits.C:69
 DrawHits.C:70
 DrawHits.C:71
 DrawHits.C:72
 DrawHits.C:73
 DrawHits.C:74
 DrawHits.C:75
 DrawHits.C:76
 DrawHits.C:77
 DrawHits.C:78
 DrawHits.C:79
 DrawHits.C:80
 DrawHits.C:81
 DrawHits.C:82
 DrawHits.C:83
 DrawHits.C:84
 DrawHits.C:85
 DrawHits.C:86
 DrawHits.C:87
 DrawHits.C:88
 DrawHits.C:89
 DrawHits.C:90
 DrawHits.C:91
 DrawHits.C:92
 DrawHits.C:93
 DrawHits.C:94
 DrawHits.C:95
 DrawHits.C:96
 DrawHits.C:97
 DrawHits.C:98
 DrawHits.C:99
 DrawHits.C:100
 DrawHits.C:101
 DrawHits.C:102
 DrawHits.C:103
 DrawHits.C:104
 DrawHits.C:105
 DrawHits.C:106
 DrawHits.C:107
 DrawHits.C:108
 DrawHits.C:109
 DrawHits.C:110
 DrawHits.C:111
 DrawHits.C:112
 DrawHits.C:113
 DrawHits.C:114
 DrawHits.C:115
 DrawHits.C:116
 DrawHits.C:117
 DrawHits.C:118
 DrawHits.C:119
 DrawHits.C:120
 DrawHits.C:121
 DrawHits.C:122
 DrawHits.C:123
 DrawHits.C:124
 DrawHits.C:125
 DrawHits.C:126
 DrawHits.C:127
 DrawHits.C:128
 DrawHits.C:129
 DrawHits.C:130
 DrawHits.C:131
 DrawHits.C:132
 DrawHits.C:133
 DrawHits.C:134
 DrawHits.C:135
 DrawHits.C:136
 DrawHits.C:137
 DrawHits.C:138
 DrawHits.C:139
 DrawHits.C:140
 DrawHits.C:141
 DrawHits.C:142
 DrawHits.C:143
 DrawHits.C:144
 DrawHits.C:145
 DrawHits.C:146
 DrawHits.C:147
 DrawHits.C:148
 DrawHits.C:149
 DrawHits.C:150
 DrawHits.C:151
 DrawHits.C:152
 DrawHits.C:153
 DrawHits.C:154
 DrawHits.C:155
 DrawHits.C:156
 DrawHits.C:157
 DrawHits.C:158
 DrawHits.C:159
 DrawHits.C:160
 DrawHits.C:161
 DrawHits.C:162
 DrawHits.C:163
 DrawHits.C:164
 DrawHits.C:165
 DrawHits.C:166
 DrawHits.C:167
 DrawHits.C:168
 DrawHits.C:169
 DrawHits.C:170
 DrawHits.C:171
 DrawHits.C:172
 DrawHits.C:173
 DrawHits.C:174
 DrawHits.C:175
 DrawHits.C:176
 DrawHits.C:177
 DrawHits.C:178
 DrawHits.C:179
 DrawHits.C:180
 DrawHits.C:181
 DrawHits.C:182
 DrawHits.C:183
 DrawHits.C:184
 DrawHits.C:185
 DrawHits.C:186
 DrawHits.C:187
 DrawHits.C:188
 DrawHits.C:189
 DrawHits.C:190
 DrawHits.C:191
 DrawHits.C:192
 DrawHits.C:193
 DrawHits.C:194
 DrawHits.C:195
 DrawHits.C:196
 DrawHits.C:197
 DrawHits.C:198
 DrawHits.C:199
 DrawHits.C:200
 DrawHits.C:201
 DrawHits.C:202
 DrawHits.C:203
 DrawHits.C:204
 DrawHits.C:205
 DrawHits.C:206
 DrawHits.C:207
 DrawHits.C:208
 DrawHits.C:209
 DrawHits.C:210
 DrawHits.C:211
 DrawHits.C:212
 DrawHits.C:213
 DrawHits.C:214
 DrawHits.C:215
 DrawHits.C:216
 DrawHits.C:217
 DrawHits.C:218
 DrawHits.C:219
 DrawHits.C:220
 DrawHits.C:221
 DrawHits.C:222
 DrawHits.C:223
 DrawHits.C:224
 DrawHits.C:225
 DrawHits.C:226
 DrawHits.C:227
 DrawHits.C:228
 DrawHits.C:229
 DrawHits.C:230
 DrawHits.C:231
 DrawHits.C:232
 DrawHits.C:233
 DrawHits.C:234
 DrawHits.C:235
 DrawHits.C:236
 DrawHits.C:237
 DrawHits.C:238
 DrawHits.C:239
 DrawHits.C:240
 DrawHits.C:241
 DrawHits.C:242
 DrawHits.C:243
 DrawHits.C:244
 DrawHits.C:245
 DrawHits.C:246
 DrawHits.C:247
 DrawHits.C:248
 DrawHits.C:249
 DrawHits.C:250
 DrawHits.C:251
 DrawHits.C:252
 DrawHits.C:253
 DrawHits.C:254
 DrawHits.C:255
 DrawHits.C:256
 DrawHits.C:257
 DrawHits.C:258
 DrawHits.C:259
 DrawHits.C:260
 DrawHits.C:261
 DrawHits.C:262
 DrawHits.C:263
 DrawHits.C:264
 DrawHits.C:265
 DrawHits.C:266
 DrawHits.C:267
 DrawHits.C:268
 DrawHits.C:269
 DrawHits.C:270
 DrawHits.C:271
 DrawHits.C:272
 DrawHits.C:273
 DrawHits.C:274
 DrawHits.C:275
 DrawHits.C:276
 DrawHits.C:277
 DrawHits.C:278
 DrawHits.C:279
 DrawHits.C:280
 DrawHits.C:281
 DrawHits.C:282
 DrawHits.C:283
 DrawHits.C:284
 DrawHits.C:285
 DrawHits.C:286
 DrawHits.C:287
 DrawHits.C:288
 DrawHits.C:289
 DrawHits.C:290
 DrawHits.C:291
 DrawHits.C:292
 DrawHits.C:293
 DrawHits.C:294
 DrawHits.C:295
 DrawHits.C:296
 DrawHits.C:297
 DrawHits.C:298
 DrawHits.C:299
 DrawHits.C:300
 DrawHits.C:301
 DrawHits.C:302
 DrawHits.C:303
 DrawHits.C:304
 DrawHits.C:305
 DrawHits.C:306
 DrawHits.C:307
 DrawHits.C:308
 DrawHits.C:309
 DrawHits.C:310
 DrawHits.C:311
 DrawHits.C:312
 DrawHits.C:313
 DrawHits.C:314
 DrawHits.C:315
 DrawHits.C:316
 DrawHits.C:317
 DrawHits.C:318
 DrawHits.C:319
 DrawHits.C:320
 DrawHits.C:321
 DrawHits.C:322
 DrawHits.C:323
 DrawHits.C:324
 DrawHits.C:325
 DrawHits.C:326
 DrawHits.C:327
 DrawHits.C:328
 DrawHits.C:329
 DrawHits.C:330
 DrawHits.C:331
 DrawHits.C:332
 DrawHits.C:333
 DrawHits.C:334
 DrawHits.C:335
 DrawHits.C:336
 DrawHits.C:337
 DrawHits.C:338
 DrawHits.C:339
 DrawHits.C:340
 DrawHits.C:341
 DrawHits.C:342
 DrawHits.C:343
 DrawHits.C:344
 DrawHits.C:345
 DrawHits.C:346
 DrawHits.C:347
 DrawHits.C:348
 DrawHits.C:349
 DrawHits.C:350
 DrawHits.C:351
 DrawHits.C:352
 DrawHits.C:353
 DrawHits.C:354
 DrawHits.C:355
 DrawHits.C:356
 DrawHits.C:357
 DrawHits.C:358
 DrawHits.C:359
 DrawHits.C:360
 DrawHits.C:361
 DrawHits.C:362
 DrawHits.C:363
 DrawHits.C:364
 DrawHits.C:365
 DrawHits.C:366
 DrawHits.C:367
 DrawHits.C:368
 DrawHits.C:369
 DrawHits.C:370
 DrawHits.C:371
 DrawHits.C:372
 DrawHits.C:373
 DrawHits.C:374
 DrawHits.C:375
 DrawHits.C:376
 DrawHits.C:377
 DrawHits.C:378
 DrawHits.C:379
 DrawHits.C:380
 DrawHits.C:381
 DrawHits.C:382
 DrawHits.C:383
 DrawHits.C:384
 DrawHits.C:385
 DrawHits.C:386
 DrawHits.C:387
 DrawHits.C:388
 DrawHits.C:389
 DrawHits.C:390
 DrawHits.C:391
 DrawHits.C:392
 DrawHits.C:393
 DrawHits.C:394
 DrawHits.C:395
 DrawHits.C:396
 DrawHits.C:397
 DrawHits.C:398
 DrawHits.C:399
 DrawHits.C:400
 DrawHits.C:401
 DrawHits.C:402
 DrawHits.C:403
 DrawHits.C:404
 DrawHits.C:405
 DrawHits.C:406
 DrawHits.C:407
 DrawHits.C:408
 DrawHits.C:409
 DrawHits.C:410
 DrawHits.C:411
 DrawHits.C:412
 DrawHits.C:413
 DrawHits.C:414
 DrawHits.C:415
 DrawHits.C:416
 DrawHits.C:417
 DrawHits.C:418
 DrawHits.C:419
 DrawHits.C:420
 DrawHits.C:421
 DrawHits.C:422
 DrawHits.C:423
 DrawHits.C:424
 DrawHits.C:425
 DrawHits.C:426
 DrawHits.C:427
 DrawHits.C:428
 DrawHits.C:429
 DrawHits.C:430
 DrawHits.C:431
 DrawHits.C:432
 DrawHits.C:433
 DrawHits.C:434
 DrawHits.C:435
 DrawHits.C:436
 DrawHits.C:437
 DrawHits.C:438
 DrawHits.C:439
 DrawHits.C:440
 DrawHits.C:441
 DrawHits.C:442
 DrawHits.C:443
 DrawHits.C:444
 DrawHits.C:445
 DrawHits.C:446
 DrawHits.C:447
 DrawHits.C:448
 DrawHits.C:449
 DrawHits.C:450
 DrawHits.C:451
 DrawHits.C:452
 DrawHits.C:453
 DrawHits.C:454
 DrawHits.C:455
 DrawHits.C:456
 DrawHits.C:457
 DrawHits.C:458
 DrawHits.C:459
 DrawHits.C:460
 DrawHits.C:461
 DrawHits.C:462
 DrawHits.C:463
 DrawHits.C:464
 DrawHits.C:465
 DrawHits.C:466
 DrawHits.C:467
 DrawHits.C:468
 DrawHits.C:469
 DrawHits.C:470
 DrawHits.C:471
 DrawHits.C:472
 DrawHits.C:473
 DrawHits.C:474
 DrawHits.C:475
 DrawHits.C:476
 DrawHits.C:477
 DrawHits.C:478
 DrawHits.C:479
 DrawHits.C:480
 DrawHits.C:481
 DrawHits.C:482
 DrawHits.C:483
 DrawHits.C:484
 DrawHits.C:485
 DrawHits.C:486
 DrawHits.C:487
 DrawHits.C:488
 DrawHits.C:489
 DrawHits.C:490
 DrawHits.C:491
 DrawHits.C:492
 DrawHits.C:493
 DrawHits.C:494
 DrawHits.C:495
 DrawHits.C:496
 DrawHits.C:497
 DrawHits.C:498
 DrawHits.C:499
 DrawHits.C:500
 DrawHits.C:501
 DrawHits.C:502
 DrawHits.C:503
 DrawHits.C:504
 DrawHits.C:505
 DrawHits.C:506
 DrawHits.C:507
 DrawHits.C:508
 DrawHits.C:509
 DrawHits.C:510
 DrawHits.C:511
 DrawHits.C:512
 DrawHits.C:513
 DrawHits.C:514
 DrawHits.C:515
 DrawHits.C:516
 DrawHits.C:517
 DrawHits.C:518
 DrawHits.C:519
 DrawHits.C:520
 DrawHits.C:521
 DrawHits.C:522
 DrawHits.C:523
 DrawHits.C:524
 DrawHits.C:525
 DrawHits.C:526
 DrawHits.C:527
 DrawHits.C:528
 DrawHits.C:529
 DrawHits.C:530
 DrawHits.C:531
 DrawHits.C:532
 DrawHits.C:533
 DrawHits.C:534
 DrawHits.C:535
 DrawHits.C:536
 DrawHits.C:537
 DrawHits.C:538
 DrawHits.C:539
 DrawHits.C:540
 DrawHits.C:541
 DrawHits.C:542
 DrawHits.C:543
 DrawHits.C:544
 DrawHits.C:545
 DrawHits.C:546
 DrawHits.C:547
 DrawHits.C:548
 DrawHits.C:549
 DrawHits.C:550
 DrawHits.C:551
 DrawHits.C:552
 DrawHits.C:553
 DrawHits.C:554
 DrawHits.C:555
 DrawHits.C:556
 DrawHits.C:557
 DrawHits.C:558
 DrawHits.C:559
 DrawHits.C:560
 DrawHits.C:561
 DrawHits.C:562
 DrawHits.C:563
 DrawHits.C:564
 DrawHits.C:565
 DrawHits.C:566
 DrawHits.C:567
 DrawHits.C:568
 DrawHits.C:569
 DrawHits.C:570
 DrawHits.C:571
 DrawHits.C:572
 DrawHits.C:573
 DrawHits.C:574
 DrawHits.C:575
 DrawHits.C:576
 DrawHits.C:577
 DrawHits.C:578
 DrawHits.C:579
 DrawHits.C:580
 DrawHits.C:581
 DrawHits.C:582
 DrawHits.C:583
 DrawHits.C:584
 DrawHits.C:585
 DrawHits.C:586
 DrawHits.C:587
 DrawHits.C:588
 DrawHits.C:589
 DrawHits.C:590
 DrawHits.C:591
 DrawHits.C:592
 DrawHits.C:593
 DrawHits.C:594
 DrawHits.C:595
 DrawHits.C:596
 DrawHits.C:597
 DrawHits.C:598
 DrawHits.C:599
 DrawHits.C:600
 DrawHits.C:601
 DrawHits.C:602
 DrawHits.C:603
 DrawHits.C:604
 DrawHits.C:605
 DrawHits.C:606
 DrawHits.C:607
 DrawHits.C:608
 DrawHits.C:609
 DrawHits.C:610
 DrawHits.C:611
 DrawHits.C:612
 DrawHits.C:613
 DrawHits.C:614
 DrawHits.C:615
 DrawHits.C:616
 DrawHits.C:617
 DrawHits.C:618
 DrawHits.C:619
 DrawHits.C:620
 DrawHits.C:621
 DrawHits.C:622
 DrawHits.C:623
 DrawHits.C:624
 DrawHits.C:625
 DrawHits.C:626
 DrawHits.C:627
 DrawHits.C:628
 DrawHits.C:629
 DrawHits.C:630
 DrawHits.C:631
 DrawHits.C:632
 DrawHits.C:633
 DrawHits.C:634
 DrawHits.C:635
 DrawHits.C:636
 DrawHits.C:637
 DrawHits.C:638
 DrawHits.C:639
 DrawHits.C:640
 DrawHits.C:641
 DrawHits.C:642
 DrawHits.C:643
 DrawHits.C:644
 DrawHits.C:645
 DrawHits.C:646
 DrawHits.C:647
 DrawHits.C:648
 DrawHits.C:649
 DrawHits.C:650
 DrawHits.C:651
 DrawHits.C:652
 DrawHits.C:653
 DrawHits.C:654
 DrawHits.C:655
 DrawHits.C:656
 DrawHits.C:657
 DrawHits.C:658
 DrawHits.C:659
 DrawHits.C:660
 DrawHits.C:661
 DrawHits.C:662
 DrawHits.C:663
 DrawHits.C:664
 DrawHits.C:665
 DrawHits.C:666
 DrawHits.C:667
 DrawHits.C:668
 DrawHits.C:669
 DrawHits.C:670
 DrawHits.C:671
 DrawHits.C:672
 DrawHits.C:673
 DrawHits.C:674
 DrawHits.C:675
 DrawHits.C:676
 DrawHits.C:677
 DrawHits.C:678
 DrawHits.C:679
 DrawHits.C:680
 DrawHits.C:681
 DrawHits.C:682
 DrawHits.C:683
 DrawHits.C:684
 DrawHits.C:685
 DrawHits.C:686
 DrawHits.C:687
 DrawHits.C:688
 DrawHits.C:689
 DrawHits.C:690
 DrawHits.C:691
 DrawHits.C:692
 DrawHits.C:693
 DrawHits.C:694
 DrawHits.C:695
 DrawHits.C:696
 DrawHits.C:697
 DrawHits.C:698
 DrawHits.C:699
 DrawHits.C:700
 DrawHits.C:701
 DrawHits.C:702
 DrawHits.C:703
 DrawHits.C:704
 DrawHits.C:705
 DrawHits.C:706
 DrawHits.C:707
 DrawHits.C:708
 DrawHits.C:709
 DrawHits.C:710
 DrawHits.C:711
 DrawHits.C:712
 DrawHits.C:713
 DrawHits.C:714
 DrawHits.C:715
 DrawHits.C:716
 DrawHits.C:717
 DrawHits.C:718
 DrawHits.C:719
 DrawHits.C:720
 DrawHits.C:721
 DrawHits.C:722
 DrawHits.C:723
 DrawHits.C:724
 DrawHits.C:725
 DrawHits.C:726
 DrawHits.C:727
 DrawHits.C:728
 DrawHits.C:729
 DrawHits.C:730
 DrawHits.C:731
 DrawHits.C:732
 DrawHits.C:733
 DrawHits.C:734
 DrawHits.C:735
 DrawHits.C:736
 DrawHits.C:737
 DrawHits.C:738
 DrawHits.C:739
 DrawHits.C:740
 DrawHits.C:741
 DrawHits.C:742
 DrawHits.C:743
 DrawHits.C:744
 DrawHits.C:745
 DrawHits.C:746
 DrawHits.C:747
 DrawHits.C:748
 DrawHits.C:749
 DrawHits.C:750
 DrawHits.C:751
 DrawHits.C:752
 DrawHits.C:753
 DrawHits.C:754
 DrawHits.C:755
 DrawHits.C:756
 DrawHits.C:757
 DrawHits.C:758