ROOT logo
/**************************************************************************
* Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
*                                                                        *
* Author: The ALICE Off-line Project.                                    *
* Contributors are mentioned in the code where appropriate.              *
*                                                                        *
* Permission to use, copy, modify and distribute this software and its   *
* documentation strictly for non-commercial purposes is hereby granted   *
* without fee, provided that the above copyright notice appears in all   *
* copies and that both the copyright notice and this permission notice   *
* appear in the supporting documentation. The authors make no claims     *
* about the suitability of this software for any purpose. It is          *
* provided "as is" without express or implied warranty.                  *
**************************************************************************/

// $Id$

/// \ingroup macros
/// \file MUONRecoCheck.C
/// \brief Utility macro to check the muon reconstruction. 
///
/// Reconstructed tracks are compared to reference tracks. The reference tracks 
/// are built from AliTrackReference for the hit in chamber (0..9) and from 
/// kinematics (TreeK) for the vertex parameters.  
///
/// \author Jean-Pierre Cussonneau, Philippe Pillot, Subatech  

// ROOT includes
#include <Riostream.h>
#include "TMath.h"
#include "TObjArray.h"
#include "TH1.h"
#include "TH2.h"
#include "TH3.h"
#include "TGraphErrors.h"
#include "TGraphAsymmErrors.h"
#include "TF1.h"
#include "TFile.h"
#include "TCanvas.h"
#include "TLegend.h"
#include "TGeoManager.h"

// STEER includes
#include "AliCDBManager.h"
#include "AliGeomManager.h"
#include "AliLog.h"

// MUON includes
#include "AliMUONCDB.h"
#include "AliMUONConstants.h"
#include "AliMUONTrack.h"
#include "AliMUONRecoCheck.h"
#include "AliMUONTrackParam.h"
#include "AliMUONRecoParam.h"
#include "AliMUONVTrackStore.h"
#include "AliMUONVCluster.h"
#include "AliMUONTrackExtrap.h"
#include "AliMUONESDInterface.h"
#include "AliMUONVTriggerTrackStore.h"
#include "AliMUONTriggerTrack.h"
#include "AliMpDEIterator.h"

Double_t langaufun(Double_t *x, Double_t *par);
void     FitGausResVsMom(TH2* h, Int_t nBins, const Double_t mean0, const Double_t sigma0, const char* fitting, TGraphAsymmErrors* gMean, TGraphAsymmErrors* gSigma);
void     FitPDCAVsMom(TH2* h, Int_t nBins, const char* fitting, TGraphAsymmErrors* gMean, TGraphAsymmErrors* gSigma);
void     FillResidual(TH1* h, Int_t i, Double_t& sigma, TGraphErrors* gMean, TGraphErrors* gSigma, Bool_t correctForSystematics, Bool_t fitResiduals);
TCanvas* DrawVsAng(const char* name, const char* title, TH1* h1, TH2* h2);
TCanvas* DrawVsPos(const char* name, const char* title, TH2* h1, TH2* h2, TH2* h3);
TCanvas* DrawResMomVsMom(const char* name, const char* title, TH2* h, Int_t nBins, TF1* f2 = 0x0, const char* fitting = "");
void     Zoom(TH1* h, Double_t fractionCut = 0.01);

//------------------------------------------------------------------------------------
void MUONRecoCheck(Int_t nEvent = -1, const char* pathSim="./generated/", const char* esdFileName="AliESDs.root",
		   const char* ocdbPath = "local://$ALICE_ROOT/OCDB", const Double_t pMin = 0., const Double_t pMax = 300.,
		   const Int_t pNBins = 30, Int_t absorberRegion = -1, Bool_t fitResiduals = kTRUE)
{
  /// Associate the reconstructed tracks with the simulated ones and check the quality of the reconstruction
  /// (tracking/trigger efficiency; momentum, slope,... resolutions at first cluster and at vertex; cluster resolution).
  /// - You can choose the momentum range and number of bins used to study the track resolution versus momentum.
  /// - You can limit the calculation of track resolution at vertex to the tracks crossing the absorber in a given region
  /// with the flag "absorberRegion": -1=all, 1=[2,3]deg, 2=[3,10]deg.
  
  Double_t aAbsLimits[2];
  if (absorberRegion > -1) {
    if (absorberRegion == 1) {
      aAbsLimits[0] = 2.;
      aAbsLimits[1] = 3.;
    } else if (absorberRegion == 2) {
      aAbsLimits[0] = 3.;
      aAbsLimits[1] = 10.;
    } else {
      cout<<"Unknown absorber region. Valid choices are: -1=all, 1=[2,3]deg, 2=[3,10]deg"<<endl;
      return;
    }
  } else {
    aAbsLimits[0] = 0.;
    aAbsLimits[1] = 90.;
  }
  
  cout<<endl<<"Momentum range for track resolution studies: "<<pNBins<<" bins in ["<<pMin<<","<<pMax<<"] GeV/c"<<endl;
  if (pMin >= pMax || pNBins <= 0) {
    cout<<"--> incorrect momentum range"<<endl<<endl;
    return;
  } else cout<<endl;
  
  AliLog::SetClassDebugLevel("AliMCEvent",-1);
  
  // ###################################### initialize ###################################### //
  AliMUONRecoCheck rc(esdFileName, pathSim);
  
  // load necessary data from OCDB
  AliCDBManager::Instance()->SetDefaultStorage(ocdbPath);
  AliCDBManager::Instance()->SetSpecificStorage("GRP/GRP/Data", Form("local://%s",gSystem->pwd()));
  AliCDBManager::Instance()->SetRun(rc.GetRunNumber());
  if (!AliMUONCDB::LoadField()) return;
  AliMUONTrackExtrap::SetField();
  if (!AliMUONCDB::LoadMapping()) return;
//  AliGeomManager::LoadGeometry();
  AliGeomManager::LoadGeometry("geometry.root");
  if (!AliGeomManager::GetGeometry()) return;
  AliMUONRecoParam* recoParam = AliMUONCDB::LoadRecoParam();
  if (!recoParam) return;
  AliMUONESDInterface::ResetTracker(recoParam);
  
  // get sigma cut from recoParam to associate clusters with TrackRefs in case the label are not used
  Double_t sigmaCut = (recoParam->ImproveTracks()) ? recoParam->GetSigmaCutForImprovement() : recoParam->GetSigmaCutForTracking();
  // compute the mask of requested stations from recoParam
  UInt_t requestedStationMask = 0;
  for (Int_t i = 0; i < 5; i++) if (recoParam->RequestStation(i)) requestedStationMask |= ( 1 << i );
  // get from recoParam whether a track need 2 chambers hit in the same station (4 or 5) or not to be reconstructible
  Bool_t request2ChInSameSt45 = !recoParam->MakeMoreTrackCandidates();
  
  // ###################################### define histograms ###################################### //
  // File for histograms and histogram booking
  TFile *histoFile = new TFile("MUONRecoCheck.root", "RECREATE");
  
  TH1F *hReconstructible = new TH1F("hReconstructible"," Nb of reconstructible tracks / evt",15,-0.5,14.5);
  TH1F *hReco = new TH1F("hReco"," Nb of reconstructed tracks / evt",15,-0.5,14.5);
  TH1F *hNClusterComp = new TH1F("hNClusterComp"," Nb of compatible clusters / track ",15,-0.5,14.5);
  TH1F *hTrackRefID = new TH1F("hTrackRefID"," track reference ID ",100,-0.5,99.5);
  TH1F *hTriggerable = new TH1F("hTriggerable"," Nb of triggerable tracks / evt",15,-0.5,14.5);
  TH1F *hTriggered = new TH1F("hTriggered"," Nb of triggered tracks / evt",15,-0.5,14.5);
  
  // momentum resolution at vertex
  histoFile->mkdir("momentumAtVertex","momentumAtVertex");
  histoFile->cd("momentumAtVertex");
  
  const Double_t pEdges[2] = {pMin, pMax};
  const Int_t deltaPAtVtxNBins = 250;
  Double_t deltaPAtVtxEdges[2];
  if (pMax < 50.) { deltaPAtVtxEdges[0] = -20.; deltaPAtVtxEdges[1] = 5.; }
  else { deltaPAtVtxEdges[0] = -50.; deltaPAtVtxEdges[1] = 50.; }
  
  TH1F *hResMomVertex = new TH1F("hResMomVertex"," delta P at vertex;#Delta_{p} (GeV/c)",deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
  
  TH2D *hResMomVertexVsMom = new TH2D("hResMomVertexVsMom","#Delta_{p} at vertex versus p;p (GeV/c);#Delta_{p} (GeV/c)",2*pNBins,pEdges[0],pEdges[1],deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
  TH2D *hResMomVertexVsMom_2_3_Deg = new TH2D("hResMomVertexVsMom_2_3_Deg","#Delta_{p} at vertex versus p for tracks between 2 and 3 degrees at absorber end;p (GeV/c);#Delta_{p} (GeV/c)",2*pNBins,pEdges[0],pEdges[1],deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
  TH2D *hResMomVertexVsMom_3_10_Deg = new TH2D("hResMomVertexVsMom_3_10_Deg","#Delta_{p} at vertex versus p for tracks between 3 and 10 degrees at absorber end;p (GeV/c);#Delta_{p} (GeV/c)",2*pNBins,pEdges[0],pEdges[1],deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
  TH2D *hResMomVertexVsMom_0_2_DegMC = new TH2D("hResMomVertexVsMom_0_2_DegMC","#Delta_{p} at vertex versus p for tracks with MC angle below 2 degrees;p (GeV/c);#Delta_{p} (GeV/c)",2*pNBins,pEdges[0],pEdges[1],deltaPAtVtxNBins/10,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
  
  TH2D *hResMomVertexVsPosAbsEnd_0_2_DegMC = new TH2D("hResMomVertexVsPosAbsEnd_0_2_DegMC","#Delta_{p} at vertex versus track position at absorber end for tracks with MC angle < 2 degrees;position (cm);#Delta_{p} (GeV/c)",1000,0.,100.,deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
  TH2D *hResMomVertexVsPosAbsEnd_2_3_DegMC = new TH2D("hResMomVertexVsPosAbsEnd_2_3_DegMC","#Delta_{p} at vertex versus track position at absorber end for tracks with MC angle in [2,3[ degrees;position (cm);#Delta_{p} (GeV/c)",1000,0.,100.,deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
  TH2D *hResMomVertexVsPosAbsEnd_3_10_DegMC = new TH2D("hResMomVertexVsPosAbsEnd_3_10_DegMC","#Delta_{p} at vertex versus track position at absorber end for tracks with MC angle in [3,10[ degrees;position (cm);#Delta_{p} (GeV/c)",1000,0.,100.,deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
  
  TH2D *hResMomVertexVsAngle = new TH2D("hResMomVertexVsAngle","#Delta_{p} at vertex versus track position at absorber end converted to degrees;angle (Deg);#Delta_{p} (GeV/c)",10,0.,10.,deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
  TH2D *hResMomVertexVsMCAngle = new TH2D("hResMomVertexVsMCAngle","#Delta_{p} at vertex versus MC angle;MC angle (Deg);#Delta_{p} (GeV/c)",10,0.,10.,deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
  TH3D *hResMomVertexVsAngleVsMom = new TH3D("hResMomVertexVsAngleVsMom","#Delta_{p} at vertex versus track position at absorber end converted to degrees versus momentum;p (GeV/c);angle (Deg);#Delta_{p} (GeV/c)",2*pNBins,pEdges[0],pEdges[1],100,0.,10.,deltaPAtVtxNBins,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1]);
  
  TGraphAsymmErrors* gMeanResMomVertexVsMom = new TGraphAsymmErrors(pNBins);
  gMeanResMomVertexVsMom->SetName("gMeanResMomVertexVsMom");
  gMeanResMomVertexVsMom->SetTitle("<#Delta_{p}> at vertex versus p;p (GeV/c);<#Delta_{p}> (GeV/c)");
  TGraphAsymmErrors* gMostProbResMomVertexVsMom = new TGraphAsymmErrors(pNBins);
  gMostProbResMomVertexVsMom->SetName("gMostProbResMomVertexVsMom");
  gMostProbResMomVertexVsMom->SetTitle("Most probable #Delta_{p} at vertex versus p;p (GeV/c);Most prob. #Delta_{p} (GeV/c)");
  TGraphAsymmErrors* gSigmaResMomVertexVsMom = new TGraphAsymmErrors(pNBins);
  gSigmaResMomVertexVsMom->SetName("gSigmaResMomVertexVsMom");
  gSigmaResMomVertexVsMom->SetTitle("#sigma_{p}/p at vertex versus p;p (GeV/c);#sigma_{p}/p (%)");
  
  // transverse momentum resolution at vertex
  TH2D *hResPtVertexVsPt = new TH2D("hResPtVertexVsPt","#Delta_{p_{t}} at vertex versus p_{t};p_{t} (GeV/c);#Delta_{p_{t}} (GeV/c)",2*pNBins,pEdges[0]/10.,pEdges[1]/10.,deltaPAtVtxNBins,deltaPAtVtxEdges[0]/10.,deltaPAtVtxEdges[1]/10.);
  
  // momentum resolution at first cluster
  histoFile->mkdir("momentumAtFirstCluster","momentumAtFirstCluster");
  histoFile->cd("momentumAtFirstCluster");
  
  const Int_t deltaPAtFirstClNBins = 500;
  Double_t deltaPAtFirstClEdges[2];
  if (pMax < 25.) { deltaPAtFirstClEdges[0] = -5.; deltaPAtFirstClEdges[1] = 5.; }
  else if (pMax < 50.) { deltaPAtFirstClEdges[0] = -10.; deltaPAtFirstClEdges[1] = 10.; }
  else { deltaPAtFirstClEdges[0] = -25.; deltaPAtFirstClEdges[1] = 25.; }
  
  TH1F *hResMomFirstCluster = new TH1F("hResMomFirstCluster"," delta P at first cluster;#Delta_{p} (GeV/c)",deltaPAtFirstClNBins,deltaPAtFirstClEdges[0],deltaPAtFirstClEdges[1]);
  TH2D *hResMomFirstClusterVsMom = new TH2D("hResMomFirstClusterVsMom","#Delta_{p} at first cluster versus p;p (GeV/c);#Delta_{p} (GeV/c)",2*pNBins,pEdges[0],pEdges[1],deltaPAtFirstClNBins,deltaPAtFirstClEdges[0],deltaPAtFirstClEdges[1]);
  
  TGraphAsymmErrors* gMeanResMomFirstClusterVsMom = new TGraphAsymmErrors(pNBins);
  gMeanResMomFirstClusterVsMom->SetName("gMeanResMomFirstClusterVsMom");
  gMeanResMomFirstClusterVsMom->SetTitle("<#Delta_{p}> at first cluster versus p;p (GeV/c);<#Delta_{p}> (GeV/c)");
  TGraphAsymmErrors* gSigmaResMomFirstClusterVsMom = new TGraphAsymmErrors(pNBins);
  gSigmaResMomFirstClusterVsMom->SetName("gSigmaResMomFirstClusterVsMom");
  gSigmaResMomFirstClusterVsMom->SetTitle("#sigma_{p}/p at first cluster versus p;p (GeV/c);#sigma_{p}/p (%)");
  
  // transverse momentum resolution at vertex
  TH2D *hResPtFirstClusterVsPt = new TH2D("hResPtFirstClusterVsPt","#Delta_{p_{t}} at first cluster versus p_{t};p_{t} (GeV/c);#Delta_{p_{t}} (GeV/c)",2*pNBins,pEdges[0]/10.,pEdges[1]/10.,deltaPAtFirstClNBins,deltaPAtFirstClEdges[0]/10.,deltaPAtFirstClEdges[1]/10.);
  
  // angular resolution at vertex
  histoFile->mkdir("slopesAtVertex","slopesAtVertex");
  histoFile->cd("slopesAtVertex");
  
  const Int_t deltaSlopeAtVtxNBins = 500;
  const Double_t deltaSlopeAtVtxEdges[2] = {-0.05, 0.05};
  
  TH1F *hResSlopeXVertex = new TH1F("hResSlopeXVertex","#Delta_{slope_{X}} at vertex;#Delta_{slope_{X}}", deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
  TH1F *hResSlopeYVertex = new TH1F("hResSlopeYVertex","#Delta_{slope_{Y}} at vertex;#Delta_{slope_{Y}}", deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
  TH2D *hResSlopeXVertexVsMom = new TH2D("hResSlopeXVertexVsMom","#Delta_{slope_{X}} at vertex versus p;p (GeV/c);#Delta_{slope_{X}}",2*pNBins,pEdges[0],pEdges[1], deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
  TH2D *hResSlopeYVertexVsMom = new TH2D("hResSlopeYVertexVsMom","#Delta_{slope_{Y}} at vertex versus p;p (GeV/c);#Delta_{slope_{Y}}",2*pNBins,pEdges[0],pEdges[1], deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
  
  TH2D *hResSlopeXVertexVsPosAbsEnd_0_2_DegMC = new TH2D("hResSlopeXVertexVsPosAbsEnd_0_2_DegMC","#Delta_{slope_{X}} at vertex versus track position at absorber end for tracks with MC angle < 2 degrees;position (cm);#Delta_{slope_{X}}",1000,0.,100.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
  TH2D *hResSlopeYVertexVsPosAbsEnd_0_2_DegMC = new TH2D("hResSlopeYVertexVsPosAbsEnd_0_2_DegMC","#Delta_{slope_{Y}} at vertex versus track position at absorber end for tracks with MC angle < 2 degrees;position (cm);#Delta_{slope_{Y}}",1000,0.,100.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
  TH2D *hResSlopeXVertexVsPosAbsEnd_2_3_DegMC = new TH2D("hResSlopeXVertexVsPosAbsEnd_2_3_DegMC","#Delta_{slope_{X}} at vertex versus track position at absorber end for tracks with MC angle in [2,3[ degrees;position (cm);#Delta_{slope_{X}}",1000,0.,100.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
  TH2D *hResSlopeYVertexVsPosAbsEnd_2_3_DegMC = new TH2D("hResSlopeYVertexVsPosAbsEnd_2_3_DegMC","#Delta_{slope_{Y}} at vertex versus track position at absorber end for tracks with MC angle in [2,3[ degrees;position (cm);#Delta_{slope_{Y}}",1000,0.,100.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
  TH2D *hResSlopeXVertexVsPosAbsEnd_3_10_DegMC = new TH2D("hResSlopeXVertexVsPosAbsEnd_3_10_DegMC","#Delta_{slope_{X}} at vertex versus track position at absorber end for tracks with MC angle in [3,10[ degrees;position (cm);#Delta_{slope_{X}}",1000,0.,100.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
  TH2D *hResSlopeYVertexVsPosAbsEnd_3_10_DegMC = new TH2D("hResSlopeYVertexVsPosAbsEnd_3_10_DegMC","#Delta_{slope_{Y}} at vertex versus track position at absorber end for tracks with MC angle in [3,10[ degrees;position (cm);#Delta_{slope_{Y}}",1000,0.,100.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
  
  TH2D *hResSlopeXVertexVsAngle = new TH2D("hResSlopeXVertexVsAngle","#Delta_{slope_{X}} at vertex versus track position at absorber end converted to degrees;angle (Deg);#Delta_{slope_{X}}",10,0.,10.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
  TH2D *hResSlopeYVertexVsAngle = new TH2D("hResSlopeYVertexVsAngle","#Delta_{slope_{Y}} at vertex versus track position at absorber end converted to degrees;angle (Deg);#Delta_{slope_{Y}}",10,0.,10.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
  TH2D *hResSlopeXVertexVsMCAngle = new TH2D("hResSlopeXVertexVsMCAngle","#Delta_{slope_{X}} at vertex versus MC angle;MC angle (Deg);#Delta_{slope_{X}}",10,0.,10.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
  TH2D *hResSlopeYVertexVsMCAngle = new TH2D("hResSlopeYVertexVsMCAngle","#Delta_{slope_{Y}} at vertex versus MC angle;MC angle (Deg);#Delta_{slope_{Y}}",10,0.,10.,deltaSlopeAtVtxNBins, deltaSlopeAtVtxEdges[0], deltaSlopeAtVtxEdges[1]);
  
  TGraphAsymmErrors* gMeanResSlopeXVertexVsMom = new TGraphAsymmErrors(pNBins);
  gMeanResSlopeXVertexVsMom->SetName("gMeanResSlopeXVertexVsMom");
  gMeanResSlopeXVertexVsMom->SetTitle("<#Delta_{slope_{X}}> at vertex versus p;p (GeV/c);<#Delta_{slope_{X}}>");
  TGraphAsymmErrors* gSigmaResSlopeXVertexVsMom = new TGraphAsymmErrors(pNBins);
  gSigmaResSlopeXVertexVsMom->SetName("gSigmaResSlopeXVertexVsMom");
  gSigmaResSlopeXVertexVsMom->SetTitle("#sigma_{slope_{X}} at vertex versus p;p (GeV/c);#sigma_{slope_{X}}");
  TGraphAsymmErrors* gMeanResSlopeYVertexVsMom = new TGraphAsymmErrors(pNBins);
  gMeanResSlopeYVertexVsMom->SetName("gMeanResSlopeYVertexVsMom");
  gMeanResSlopeYVertexVsMom->SetTitle("<#Delta_{slope_{Y}}> at vertex versus p;p (GeV/c);<#Delta_{slope_{Y}}>");
  TGraphAsymmErrors* gSigmaResSlopeYVertexVsMom = new TGraphAsymmErrors(pNBins);
  gSigmaResSlopeYVertexVsMom->SetName("gSigmaResSlopeYVertexVsMom");
  gSigmaResSlopeYVertexVsMom->SetTitle("#sigma_{slope_{Y}} at vertex versus p;p (GeV/c);#sigma_{slope_{Y}}");
  
  // angular resolution at first cluster
  histoFile->mkdir("slopesAtFirstCluster","slopesAtFirstCluster");
  histoFile->cd("slopesAtFirstCluster");
  
  const Int_t deltaSlopeAtFirstClNBins = 500;
  const Double_t deltaSlopeAtFirstClEdges[2] = {-0.01, 0.01};
  
  TH1F *hResSlopeXFirstCluster = new TH1F("hResSlopeXFirstCluster","#Delta_{slope_{X}} at first cluster;#Delta_{slope_{X}}", deltaSlopeAtFirstClNBins, deltaSlopeAtFirstClEdges[0], deltaSlopeAtFirstClEdges[1]);
  TH2D *hResSlopeXFirstClusterVsMom = new TH2D("hResSlopeXFirstClusterVsMom","#Delta_{slope_{X}} at first cluster versus p;p (GeV/c);#Delta_{slope_{X}}",2*pNBins,pEdges[0],pEdges[1], deltaSlopeAtFirstClNBins, deltaSlopeAtFirstClEdges[0], deltaSlopeAtFirstClEdges[1]);
  TH1F *hResSlopeYFirstCluster = new TH1F("hResSlopeYFirstCluster","#Delta_{slope_{Y}} at first cluster;#Delta_{slope_{Y}}", deltaSlopeAtFirstClNBins, deltaSlopeAtFirstClEdges[0], deltaSlopeAtFirstClEdges[1]);
  TH2D *hResSlopeYFirstClusterVsMom = new TH2D("hResSlopeYFirstClusterVsMom","#Delta_{slope_{Y}} at first cluster versus p;p (GeV/c);#Delta_{slope_{Y}}",2*pNBins,pEdges[0],pEdges[1], deltaSlopeAtFirstClNBins, deltaSlopeAtFirstClEdges[0], deltaSlopeAtFirstClEdges[1]);
  
  TGraphAsymmErrors* gMeanResSlopeXFirstClusterVsMom = new TGraphAsymmErrors(pNBins);
  gMeanResSlopeXFirstClusterVsMom->SetName("gMeanResSlopeXFirstClusterVsMom");
  gMeanResSlopeXFirstClusterVsMom->SetTitle("<#Delta_{slope_{X}}> at first cluster versus p;p (GeV/c);<#Delta_{slope_{X}}>");
  TGraphAsymmErrors* gSigmaResSlopeXFirstClusterVsMom = new TGraphAsymmErrors(pNBins);
  gSigmaResSlopeXFirstClusterVsMom->SetName("gSigmaResSlopeXFirstClusterVsMom");
  gSigmaResSlopeXFirstClusterVsMom->SetTitle("#sigma_{slope_{X}} at first cluster versus p;p (GeV/c);#sigma_{slope_{X}}");
  TGraphAsymmErrors* gMeanResSlopeYFirstClusterVsMom = new TGraphAsymmErrors(pNBins);
  gMeanResSlopeYFirstClusterVsMom->SetName("gMeanResSlopeYFirstClusterVsMom");
  gMeanResSlopeYFirstClusterVsMom->SetTitle("<#Delta_{slope_{Y}}> at first cluster versus p;p (GeV/c);<#Delta_{slope_{Y}}>");
  TGraphAsymmErrors* gSigmaResSlopeYFirstClusterVsMom = new TGraphAsymmErrors(pNBins);
  gSigmaResSlopeYFirstClusterVsMom->SetName("gSigmaResSlopeYFirstClusterVsMom");
  gSigmaResSlopeYFirstClusterVsMom->SetTitle("#sigma_{slope_{Y}} at first cluster versus p;p (GeV/c);#sigma_{slope_{Y}}");
  
  // DCA resolution and MCS angular dispersion
  histoFile->mkdir("DCA","DCA");
  histoFile->cd("DCA");
  
  const Int_t deltaPDCANBins = 500;
  const Double_t deltaPDCAEdges[2] = {0., 1000.};
  const Double_t deltaPMCSAngEdges[2] = {-1.5, 1.5};
  
  TH1F *hPDCA = new TH1F("hPDCA","p #times DCA at vertex;p #times DCA (GeV #times cm)", deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
  TH2D *hPDCAVsMom_2_3_Deg = new TH2D("hPDCAVsMom_2_3_Deg","p #times DCA versus p for tracks within [2,3[ degrees at absorber end;p (GeV/c);p #times DCA (GeV #times cm)",2*pNBins,pEdges[0],pEdges[1], deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
  TH2D *hPDCAVsMom_3_10_Deg = new TH2D("hPDCAVsMom_3_10_Deg","p #times DCA versus p for tracks within [3,10[ degrees at absorber end;p (GeV/c);p #times DCA (GeV #times cm)",2*pNBins,pEdges[0],pEdges[1], deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
  TH2D *hPMCSAngVsMom_2_3_Deg = new TH2D("hPMCSAngVsMom_2_3_Deg","p #times #Delta#theta_{MCS} versus p for tracks within [2,3[ degrees at absorber end;p (GeV/c);p #times #Delta#theta_{MCS} (GeV)",2*pNBins,pEdges[0],pEdges[1], deltaPDCANBins, deltaPMCSAngEdges[0], deltaPMCSAngEdges[1]);
  TH2D *hPMCSAngVsMom_3_10_Deg = new TH2D("hPMCSAngVsMom_3_10_Deg","p #times #Delta#theta_{MCS} versus p for tracks within [2,3[ degrees at absorber end;p (GeV/c);p #times #Delta#theta_{MCS} (GeV)",2*pNBins,pEdges[0],pEdges[1], deltaPDCANBins, deltaPMCSAngEdges[0], deltaPMCSAngEdges[1]);
  
  TH2D *hPDCAVsPosAbsEnd_0_2_DegMC = new TH2D("hPDCAVsPosAbsEnd_0_2_DegMC","p #times DCA versus track position at absorber end for tracks with MC angle < 2 degrees;position (cm);p #times DCA (GeV #times cm)",1000,0.,100.,deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
  TH2D *hPDCAVsPosAbsEnd_2_3_DegMC = new TH2D("hPDCAVsPosAbsEnd_2_3_DegMC","p #times DCA}versus track position at absorber end for tracks with MC angle in [2,3[ degrees;position (cm);p #times DCA (GeV #times cm)",1000,0.,100.,deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
  TH2D *hPDCAVsPosAbsEnd_3_10_DegMC = new TH2D("hPDCAVsPosAbsEnd_3_10_DegMC","p #times DCA versus track position at absorber end for tracks with MC angle in [3,10[ degrees;position (cm);p #times DCA (GeV #times cm)",1000,0.,100.,deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
  
  TH2D *hPDCAVsAngle = new TH2D("hPDCAVsAngle","p #times DCA versus track position at absorber end converted to degrees;angle (Deg);p #times DCA (GeV #times cm)",10,0.,10.,deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
  TH2D *hPDCAVsMCAngle = new TH2D("hPDCAVsMCAngle","p #times DCA versus MC angle;MC angle (Deg);p #times DCA (GeV #times cm)",10,0.,10.,deltaPDCANBins, deltaPDCAEdges[0], deltaPDCAEdges[1]);
  
  TGraphAsymmErrors* gMeanPDCAVsMom_2_3_Deg = new TGraphAsymmErrors(pNBins);
  gMeanPDCAVsMom_2_3_Deg->SetName("gMeanPDCAVsMom_2_3_Deg");
  gMeanPDCAVsMom_2_3_Deg->SetTitle("<p #times DCA> versus p for tracks within [2,3[ degrees at absorber end;p (GeV/c);<p #times DCA> (GeV #times cm)");
  TGraphAsymmErrors* gSigmaPDCAVsMom_2_3_Deg = new TGraphAsymmErrors(pNBins);
  gSigmaPDCAVsMom_2_3_Deg->SetName("gSigmaPDCAVsMom_2_3_Deg");
  gSigmaPDCAVsMom_2_3_Deg->SetTitle("#sigma_{p #times DCA} versus p for tracks within [2,3[ degrees at absorber end;p (GeV/c);#sigma_{p #times DCA} (GeV #times cm)");
  TGraphAsymmErrors* gMeanPDCAVsMom_3_10_Deg = new TGraphAsymmErrors(pNBins);
  gMeanPDCAVsMom_3_10_Deg->SetName("gMeanPDCAVsMom_3_10_Deg");
  gMeanPDCAVsMom_3_10_Deg->SetTitle("<p #times DCA> versus p for tracks within [3,10[ degrees at absorber end;p (GeV/c);<p #times DCA> (GeV #times cm)");
  TGraphAsymmErrors* gSigmaPDCAVsMom_3_10_Deg = new TGraphAsymmErrors(pNBins);
  gSigmaPDCAVsMom_3_10_Deg->SetName("gSigmaPDCAVsMom_3_10_Deg");
  gSigmaPDCAVsMom_3_10_Deg->SetTitle("#sigma_{p #times DCA} versus p for tracks within [3,10[ degrees at absorber end;p (GeV/c);#sigma_{p #times DCA} (GeV #times cm)");
  TGraphAsymmErrors* gMeanPMCSAngVsMom_2_3_Deg = new TGraphAsymmErrors(pNBins);
  gMeanPMCSAngVsMom_2_3_Deg->SetName("gMeanPMCSAngVsMom_2_3_Deg");
  gMeanPMCSAngVsMom_2_3_Deg->SetTitle("<p #times #Delta#theta_{MCS}> versus p for tracks within [2,3[ degrees at absorber end;p (GeV/c);<p #times #Delta#theta_{MCS}> (GeV)");
  TGraphAsymmErrors* gSigmaPMCSAngVsMom_2_3_Deg = new TGraphAsymmErrors(pNBins);
  gSigmaPMCSAngVsMom_2_3_Deg->SetName("gSigmaPMCSAngVsMom_2_3_Deg");
  gSigmaPMCSAngVsMom_2_3_Deg->SetTitle("#sigma_{p #times #Delta#theta_{MCS}} versus p for tracks within [2,3[ degrees at absorber end;p (GeV/c);#sigma_{p #times #Delta#theta_{MCS}} (GeV)");
  TGraphAsymmErrors* gMeanPMCSAngVsMom_3_10_Deg = new TGraphAsymmErrors(pNBins);
  gMeanPMCSAngVsMom_3_10_Deg->SetName("gMeanPMCSAngVsMom_3_10_Deg");
  gMeanPMCSAngVsMom_3_10_Deg->SetTitle("<p #times #Delta#theta_{MCS}> versus p for tracks within [3,10[ degrees at absorber end;p (GeV/c);<p #times #Delta#theta_{MCS}> (GeV)");
  TGraphAsymmErrors* gSigmaPMCSAngVsMom_3_10_Deg = new TGraphAsymmErrors(pNBins);
  gSigmaPMCSAngVsMom_3_10_Deg->SetName("gSigmaPMCSAngVsMom_3_10_Deg");
  gSigmaPMCSAngVsMom_3_10_Deg->SetTitle("#sigma_{p #times #Delta#theta_{MCS}} versus p for tracks within [3,10[ degrees at absorber end;p (GeV/c);#sigma_{p #times #Delta#theta_{MCS}} (GeV)");
  
  // eta resolution at vertex
  histoFile->mkdir("etaAtVertex","etaAtVertex");
  histoFile->cd("etaAtVertex");
  
  const Int_t deltaEtaAtVtxNBins = 500;
  const Double_t deltaEtaAtVtxEdges[2] = {-0.5, 0.5};
  
  TH1F *hResEtaVertex = new TH1F("hResEtaVertex","#Delta_{eta} at vertex;#Delta_{eta}", deltaEtaAtVtxNBins, deltaEtaAtVtxEdges[0], deltaEtaAtVtxEdges[1]);
  TH2D *hResEtaVertexVsMom = new TH2D("hResEtaVertexVsMom","#Delta_{eta} at vertex versus p;p (GeV/c);#Delta_{eta}",2*pNBins,pEdges[0],pEdges[1], deltaEtaAtVtxNBins, deltaEtaAtVtxEdges[0], deltaEtaAtVtxEdges[1]);
  
  TH2D *hResEtaVertexVsPosAbsEnd_0_2_DegMC = new TH2D("hResEtaVertexVsPosAbsEnd_0_2_DegMC","#Delta_{eta} at vertex versus track position at absorber end for tracks with MC angle < 2 degrees;position (cm);#Delta_{eta}",1000,0.,100.,deltaEtaAtVtxNBins, deltaEtaAtVtxEdges[0], deltaEtaAtVtxEdges[1]);
  TH2D *hResEtaVertexVsPosAbsEnd_2_3_DegMC = new TH2D("hResEtaVertexVsPosAbsEnd_2_3_DegMC","#Delta_{eta} at vertex versus track position at absorber end for tracks with MC angle in [2,3[ degrees;position (cm);#Delta_{eta}",1000,0.,100.,deltaEtaAtVtxNBins, deltaEtaAtVtxEdges[0], deltaEtaAtVtxEdges[1]);
  TH2D *hResEtaVertexVsPosAbsEnd_3_10_DegMC = new TH2D("hResEtaVertexVsPosAbsEnd_3_10_DegMC","#Delta_{eta} at vertex versus track position at absorber end for tracks with MC angle in [3,10[ degrees;position (cm);#Delta_{eta}",1000,0.,100.,deltaEtaAtVtxNBins, deltaEtaAtVtxEdges[0], deltaEtaAtVtxEdges[1]);
  
  TH2D *hResEtaVertexVsAngle = new TH2D("hResEtaVertexVsAngle","#Delta_{eta} at vertex versus track position at absorber end converted to degrees;angle (Deg);#Delta_{eta}",10,0.,10.,deltaEtaAtVtxNBins, deltaEtaAtVtxEdges[0], deltaEtaAtVtxEdges[1]);
  TH2D *hResEtaVertexVsMCAngle = new TH2D("hResEtaVertexVsMCAngle","#Delta_{eta} at vertex versus MC angle;MC angle (Deg);#Delta_{eta}",10,0.,10.,deltaEtaAtVtxNBins, deltaEtaAtVtxEdges[0], deltaEtaAtVtxEdges[1]);
  
  TGraphAsymmErrors* gMeanResEtaVertexVsMom = new TGraphAsymmErrors(pNBins);
  gMeanResEtaVertexVsMom->SetName("gMeanResEtaVertexVsMom");
  gMeanResEtaVertexVsMom->SetTitle("<#Delta_{eta}> at vertex versus p;p (GeV/c);<#Delta_{eta}>");
  TGraphAsymmErrors* gSigmaResEtaVertexVsMom = new TGraphAsymmErrors(pNBins);
  gSigmaResEtaVertexVsMom->SetName("gSigmaResEtaVertexVsMom");
  gSigmaResEtaVertexVsMom->SetTitle("#sigma_{eta} at vertex versus p;p (GeV/c);#sigma_{eta}");
  
  // phi resolution at vertex
  histoFile->mkdir("phiAtVertex","phiAtVertex");
  histoFile->cd("phiAtVertex");
  
  const Int_t deltaPhiAtVtxNBins = 500;
  const Double_t deltaPhiAtVtxEdges[2] = {-0.5, 0.5};
  
  TH1F *hResPhiVertex = new TH1F("hResPhiVertex","#Delta_{phi} at vertex;#Delta_{phi}", deltaPhiAtVtxNBins, deltaPhiAtVtxEdges[0], deltaPhiAtVtxEdges[1]);
  TH2D *hResPhiVertexVsMom = new TH2D("hResPhiVertexVsMom","#Delta_{phi} at vertex versus p;p (GeV/c);#Delta_{phi}",2*pNBins,pEdges[0],pEdges[1], deltaPhiAtVtxNBins, deltaPhiAtVtxEdges[0], deltaPhiAtVtxEdges[1]);
  
  TH2D *hResPhiVertexVsPosAbsEnd_0_2_DegMC = new TH2D("hResPhiVertexVsPosAbsEnd_0_2_DegMC","#Delta_{phi} at vertex versus track position at absorber end for tracks with MC angle < 2 degrees;position (cm);#Delta_{phi}",1000,0.,100.,deltaPhiAtVtxNBins, deltaPhiAtVtxEdges[0], deltaPhiAtVtxEdges[1]);
  TH2D *hResPhiVertexVsPosAbsEnd_2_3_DegMC = new TH2D("hResPhiVertexVsPosAbsEnd_2_3_DegMC","#Delta_{phi} at vertex versus track position at absorber end for tracks with MC angle in [2,3[ degrees;position (cm);#Delta_{phi}",1000,0.,100.,deltaPhiAtVtxNBins, deltaPhiAtVtxEdges[0], deltaPhiAtVtxEdges[1]);
  TH2D *hResPhiVertexVsPosAbsEnd_3_10_DegMC = new TH2D("hResPhiVertexVsPosAbsEnd_3_10_DegMC","#Delta_{phi} at vertex versus track position at absorber end for tracks with MC angle in [3,10[ degrees;position (cm);#Delta_{phi}",1000,0.,100.,deltaPhiAtVtxNBins, deltaPhiAtVtxEdges[0], deltaPhiAtVtxEdges[1]);
  
  TH2D *hResPhiVertexVsAngle = new TH2D("hResPhiVertexVsAngle","#Delta_{phi} at vertex versus track position at absorber end converted to degrees;angle (Deg);#Delta_{phi}",10,0.,10.,deltaPhiAtVtxNBins, deltaPhiAtVtxEdges[0], deltaPhiAtVtxEdges[1]);
  TH2D *hResPhiVertexVsMCAngle = new TH2D("hResPhiVertexVsMCAngle","#Delta_{phi} at vertex versus MC angle;MC angle (Deg);#Delta_{phi}",10,0.,10.,deltaPhiAtVtxNBins, deltaPhiAtVtxEdges[0], deltaPhiAtVtxEdges[1]);
  
  TGraphAsymmErrors* gMeanResPhiVertexVsMom = new TGraphAsymmErrors(pNBins);
  gMeanResPhiVertexVsMom->SetName("gMeanResPhiVertexVsMom");
  gMeanResPhiVertexVsMom->SetTitle("<#Delta_{phi}> at vertex versus p;p (GeV/c);<#Delta_{phi}>");
  TGraphAsymmErrors* gSigmaResPhiVertexVsMom = new TGraphAsymmErrors(pNBins);
  gSigmaResPhiVertexVsMom->SetName("gSigmaResPhiVertexVsMom");
  gSigmaResPhiVertexVsMom->SetTitle("#sigma_{phi} at vertex versus p;p (GeV/c);#sigma_{phi}");
  
  // cluster resolution
  histoFile->mkdir("clusters","clusters");
  histoFile->cd("clusters");
  
  // find the highest chamber resolution and set histogram bins
  const Int_t clusterResNBins = 5000;
  Double_t maxRes[2] = {-1., -1.};
  for (Int_t i = 0; i < 10; i++) {
    if (recoParam->GetDefaultNonBendingReso(i) > maxRes[0]) maxRes[0] = recoParam->GetDefaultNonBendingReso(i);
    if (recoParam->GetDefaultBendingReso(i) > maxRes[1]) maxRes[1] = recoParam->GetDefaultBendingReso(i);
  }
  Double_t clusterResMaxX = 10.*maxRes[0];
  Double_t clusterResMaxY = 10.*maxRes[1];
  
  TH1F* hResidualXInCh[AliMUONConstants::NTrackingCh()];
  TH1F* hResidualYInCh[AliMUONConstants::NTrackingCh()];
  for (Int_t i = 0; i < AliMUONConstants::NTrackingCh(); i++) {
    hResidualXInCh[i] = new TH1F(Form("hResidualXInCh%d",i+1), Form("cluster-track residual-X distribution in chamber %d;#Delta_{X} (cm)",i+1), clusterResNBins, -clusterResMaxX, clusterResMaxX);
    hResidualYInCh[i] = new TH1F(Form("hResidualYInCh%d",i+1), Form("cluster-track residual-Y distribution in chamber %d;#Delta_{Y} (cm)",i+1), clusterResNBins, -clusterResMaxY, clusterResMaxY);
  }
  
  // fill correspondence between DEId and indices for histo (starting from 1)
  Int_t deNBins = 0;
  Int_t deIndices[1100];
  Int_t deIds[200];
  AliMpDEIterator it;
  for (Int_t i = 0; i < AliMUONConstants::NTrackingCh(); i++) {
    it.First(i);
    while (!it.IsDone()) {
      deNBins++;
      deIndices[it.CurrentDEId()] = deNBins;
      deIds[deNBins] = it.CurrentDEId();
      it.Next();
    }
  }
  TH2F* hResidualXPerDE = new TH2F("hResidualXPerDE", "cluster-track residual-X distribution per DE;DE ID;#Delta_{X} (cm)", deNBins, 0.5, deNBins+0.5, clusterResNBins, -clusterResMaxX, clusterResMaxX);
  for (Int_t i = 1; i <= deNBins; i++) hResidualXPerDE->GetXaxis()->SetBinLabel(i, Form("%d",deIds[i]));
  TH2F* hResidualYPerDE = new TH2F("hResidualYPerDE", "cluster-track residual-Y distribution per DE;DE ID;#Delta_{Y} (cm)", deNBins, 0.5, deNBins+0.5, clusterResNBins, -clusterResMaxY, clusterResMaxY);
  for (Int_t i = 1; i <= deNBins; i++) hResidualYPerDE->GetXaxis()->SetBinLabel(i, Form("%d",deIds[i]));
  
  TGraphErrors* gResidualXPerChMean = new TGraphErrors(AliMUONConstants::NTrackingCh());
  gResidualXPerChMean->SetName("gResidualXPerChMean");
  gResidualXPerChMean->SetTitle("cluster-trackRef residual-X per Ch: mean;chamber ID;<#Delta_{X}> (cm)");
  gResidualXPerChMean->SetMarkerStyle(kFullDotLarge);
  TGraphErrors* gResidualYPerChMean = new TGraphErrors(AliMUONConstants::NTrackingCh());
  gResidualYPerChMean->SetName("gResidualYPerChMean");
  gResidualYPerChMean->SetTitle("cluster-trackRef residual-Y per Ch: mean;chamber ID;<#Delta_{Y}> (cm)");
  gResidualYPerChMean->SetMarkerStyle(kFullDotLarge);
  TGraphErrors* gResidualXPerChSigma = new TGraphErrors(AliMUONConstants::NTrackingCh());
  gResidualXPerChSigma->SetName("gResidualXPerChSigma");
  gResidualXPerChSigma->SetTitle("cluster-trackRef residual-X per Ch: sigma;chamber ID;#sigma_{X} (cm)");
  gResidualXPerChSigma->SetMarkerStyle(kFullDotLarge);
  TGraphErrors* gResidualYPerChSigma = new TGraphErrors(AliMUONConstants::NTrackingCh());
  gResidualYPerChSigma->SetName("gResidualYPerChSigma");
  gResidualYPerChSigma->SetTitle("cluster-trackRef residual-Y per Ch: sigma;chamber ID;#sigma_{Y} (cm)");
  gResidualYPerChSigma->SetMarkerStyle(kFullDotLarge);
  
  TGraphErrors* gResidualXPerDEMean = new TGraphErrors(deNBins);
  gResidualXPerDEMean->SetName("gResidualXPerDEMean");
  gResidualXPerDEMean->SetTitle("cluster-trackRef residual-X per DE: mean;DE ID;<#Delta_{X}> (cm)");
  gResidualXPerDEMean->SetMarkerStyle(kFullDotLarge);
  TGraphErrors* gResidualYPerDEMean = new TGraphErrors(deNBins);
  gResidualYPerDEMean->SetName("gResidualYPerDEMean");
  gResidualYPerDEMean->SetTitle("cluster-trackRef residual-Y per dE: mean;DE ID;<#Delta_{Y}> (cm)");
  gResidualYPerDEMean->SetMarkerStyle(kFullDotLarge);
  TGraphErrors* gResidualXPerDESigma = new TGraphErrors(deNBins);
  gResidualXPerDESigma->SetName("gResidualXPerDESigma");
  gResidualXPerDESigma->SetTitle("cluster-trackRef residual-X per DE: sigma;DE ID;#sigma_{X} (cm)");
  gResidualXPerDESigma->SetMarkerStyle(kFullDotLarge);
  TGraphErrors* gResidualYPerDESigma = new TGraphErrors(deNBins);
  gResidualYPerDESigma->SetName("gResidualYPerDESigma");
  gResidualYPerDESigma->SetTitle("cluster-trackRef residual-Y per DE: sigma;DE ID;#sigma_{Y} (cm)");
  gResidualYPerDESigma->SetMarkerStyle(kFullDotLarge);
  
  histoFile->mkdir("trigger");
  histoFile->cd("trigger");
  TH1F* hResidualTrigX11 = new TH1F("hResiudalTrigX11", "Residual X11", 100, -10., 10.);
  TH1F* hResidualTrigY11 = new TH1F("hResiudalTrigY11", "Residual Y11", 100, -10., 10.);
  TH1F* hResidualTrigSlopeY = new TH1F("hResiudalTrigSlopeY", "Residual Y slope", 100, -0.1, 0.1);
  TH1F* hTriggerableMatchFailed = new TH1F("hTriggerableMatchFailed", "Triggerable multiplicity for events with no match", 15, -0.5, 14.5);
  
  // ###################################### fill histograms ###################################### //
  Int_t ievent;
  Int_t nReconstructibleTracks = 0;
  Int_t nReconstructedTracks = 0;
  Int_t nReconstructibleTracksCheck = 0;
  AliMUONTrackParam *trackParam;
  Double_t x1,y1,z1,slopex1,slopey1,pX1,pY1,pZ1,p1,pT1,eta1,phi1;
  Double_t x2,y2,z2,slopex2,slopey2,pX2,pY2,pZ2,p2,pT2,eta2,phi2;
  Double_t dPhi;
  Double_t xAbs,yAbs,dAbs,aAbs,aMCS,aMC;
  Double_t xDCA,yDCA,dca,pU;
  Double_t aMCSMoy = 0., aMCS2Moy = 0., dMCSMoy = 0., dMCS2Moy = 0., adMCSMoy = 0.;
  Int_t nMCS = 0;
  
  Int_t nevents = rc.NumberOfEvents();
  if (nevents < nEvent || nEvent < 0) nEvent = nevents;
  
  // loop over events
  for (ievent=0; ievent<nEvent; ievent++)
  {
    if ((ievent+1)%100 == 0) cout<<"\rEvent processing... "<<ievent+1<<flush;
    
    AliMUONVTrackStore* trackStore = rc.ReconstructedTracks(ievent, kFALSE);
    AliMUONVTrackStore* trackRefStore = rc.ReconstructibleTracks(ievent, requestedStationMask, request2ChInSameSt45);
    
    hReconstructible->Fill(trackRefStore->GetSize());
    hReco->Fill(trackStore->GetSize());
    
    nReconstructibleTracks += trackRefStore->GetSize();
    nReconstructedTracks += trackStore->GetSize();

    AliMUONVTriggerTrackStore* triggerTrackRefStore = rc.TriggerableTracks(ievent);
    AliMUONVTriggerTrackStore* triggerTrackStore = rc.TriggeredTracks(ievent);

    hTriggerable->Fill(triggerTrackRefStore->GetSize());
    hTriggered->Fill(triggerTrackStore->GetSize());

    // loop over trigger trackRef
    TIter nextTrig(triggerTrackRefStore->CreateIterator());
    AliMUONTriggerTrack* triggerTrackRef;
    Int_t nTriggerMatches = 0;
    while ( ( triggerTrackRef = static_cast<AliMUONTriggerTrack*>(nextTrig()) ) )
    {
      
      AliMUONTriggerTrack* triggerTrackMatched = 0x0;
      
      // loop over trackReco and look for compatible track
      TIter nextTrig2(triggerTrackStore->CreateIterator());
      AliMUONTriggerTrack* triggerTrackReco;
      while ( ( triggerTrackReco = static_cast<AliMUONTriggerTrack*>(nextTrig2()) ) )
      {
	
        // check if trackReco is compatible with trackRef
        if (triggerTrackReco->Match(*triggerTrackRef, sigmaCut)) {
          triggerTrackMatched = triggerTrackReco;
          nTriggerMatches++;
          break;
        }
      }
      
      if (triggerTrackMatched) { // tracking requirements verified, track is found
        hResidualTrigX11->Fill( triggerTrackMatched->GetX11() - triggerTrackRef->GetX11() );
        hResidualTrigY11->Fill( triggerTrackMatched->GetY11() - triggerTrackRef->GetY11() );
        hResidualTrigSlopeY->Fill( triggerTrackMatched->GetSlopeY() - triggerTrackRef->GetSlopeY() );
      }
    } // loop on trigger track ref
    
    if ( nTriggerMatches != triggerTrackStore->GetSize() )
      hTriggerableMatchFailed->Fill(triggerTrackRefStore->GetSize());
    
    // loop over trackRef
    TIter next(trackRefStore->CreateIterator());
    AliMUONTrack* trackRef;
    while ( ( trackRef = static_cast<AliMUONTrack*>(next()) ) )
    {
      
      hTrackRefID->Fill(trackRef->GetUniqueID());
      
      AliMUONTrack* trackMatched = 0x0;
      Int_t nMatchClusters = 0;
      
      // loop over trackReco and look for compatible track
      TIter next2(trackStore->CreateIterator());
      AliMUONTrack* trackReco;
      while ( ( trackReco = static_cast<AliMUONTrack*>(next2()) ) )
      {
	
	// check if trackReco is compatible with trackRef
	if (trackReco->Match(*trackRef, sigmaCut, nMatchClusters)) {
	  trackMatched = trackReco;
	  break;
	}
	
      }
      
      if (trackMatched) { // tracking requirements verified, track is found
        nReconstructibleTracksCheck++;
        hNClusterComp->Fill(nMatchClusters);
	
	// compute track position at the end of the absorber
        AliMUONTrackParam trackParamAtAbsEnd(*((AliMUONTrackParam*)trackMatched->GetTrackParamAtCluster()->First()));
	AliMUONTrackExtrap::ExtrapToZ(&trackParamAtAbsEnd, AliMUONConstants::AbsZEnd());
        xAbs = trackParamAtAbsEnd.GetNonBendingCoor();
        yAbs = trackParamAtAbsEnd.GetBendingCoor();
	dAbs = TMath::Sqrt(xAbs*xAbs + yAbs*yAbs);
	aAbs = TMath::ATan(-dAbs/AliMUONConstants::AbsZEnd()) * TMath::RadToDeg();
        pX2 = trackParamAtAbsEnd.Px();
        pY2 = trackParamAtAbsEnd.Py();
        pZ2 = trackParamAtAbsEnd.Pz();
        pT2 = TMath::Sqrt(pX2*pX2 + pY2*pY2);
	aMCS = TMath::ATan(-pT2/pZ2) * TMath::RadToDeg();
	
        trackParam = trackRef->GetTrackParamAtVertex();
        x1 = trackParam->GetNonBendingCoor();
        y1 = trackParam->GetBendingCoor();
        z1 = trackParam->GetZ();
        slopex1 = trackParam->GetNonBendingSlope();
        slopey1 = trackParam->GetBendingSlope();
        pX1 = trackParam->Px();
        pY1 = trackParam->Py();
        pZ1 = trackParam->Pz();
        p1  = trackParam->P();
        pT1 = TMath::Sqrt(pX1*pX1 + pY1*pY1);
	aMC = TMath::ATan(-pT1/pZ1) * TMath::RadToDeg();
	eta1 = TMath::Log(TMath::Tan(0.5*TMath::ATan(-pT1/pZ1)));
	phi1 = TMath::Pi()+TMath::ATan2(-pY1, -pX1);
	
	trackParam = trackMatched->GetTrackParamAtVertex();
        x2 = trackParam->GetNonBendingCoor();
        y2 = trackParam->GetBendingCoor();
        z2 = trackParam->GetZ();
        slopex2 = trackParam->GetNonBendingSlope();
        slopey2 = trackParam->GetBendingSlope();
        pX2 = trackParam->Px();
        pY2 = trackParam->Py();
        pZ2 = trackParam->Pz();
        p2  = trackParam->P();
        pT2 = TMath::Sqrt(pX2*pX2 + pY2*pY2);
	eta2 = TMath::Log(TMath::Tan(0.5*TMath::ATan(-pT2/pZ2)));
	phi2 = TMath::Pi()+TMath::ATan2(-pY2, -pX2);
        
	dPhi = phi2-phi1;
	if (dPhi < -TMath::Pi()) dPhi += 2.*TMath::Pi();
	else if (dPhi > TMath::Pi()) dPhi -= 2.*TMath::Pi();
	
        AliMUONTrackParam trackParamAtDCA(*((AliMUONTrackParam*) trackMatched->GetTrackParamAtCluster()->First()));
	pU = trackParamAtDCA.P();
	AliMUONTrackExtrap::ExtrapToVertexWithoutBranson(&trackParamAtDCA, z2);
        xDCA = trackParamAtDCA.GetNonBendingCoor();
        yDCA = trackParamAtDCA.GetBendingCoor();
	dca = TMath::Sqrt(xDCA*xDCA + yDCA*yDCA);
	
        hResMomVertex->Fill(p2-p1);
	hResSlopeXVertex->Fill(slopex2-slopex1);
	hResSlopeYVertex->Fill(slopey2-slopey1);
	hPDCA->Fill(0.5*(p2+pU)*dca);
	hResEtaVertex->Fill(eta2-eta1);
	hResPhiVertex->Fill(dPhi);
	if (aMC >= aAbsLimits[0] && aMC <= aAbsLimits[1]) {
	  hResMomVertexVsMom->Fill(p1,p2-p1);
	  hResSlopeXVertexVsMom->Fill(p1,slopex2-slopex1);
	  hResSlopeYVertexVsMom->Fill(p1,slopey2-slopey1);
	  hResEtaVertexVsMom->Fill(p1,eta2-eta1);
	  hResPhiVertexVsMom->Fill(p1,dPhi);
	  hResPtVertexVsPt->Fill(pT1,pT2-pT1);
	}
	hResMomVertexVsAngleVsMom->Fill(p1,aAbs,p2-p1);
	if (aAbs > 2. && aAbs < 3.) {
	  hResMomVertexVsMom_2_3_Deg->Fill(p1,p2-p1);
	  hPDCAVsMom_2_3_Deg->Fill(p1,0.5*(p2+pU)*dca);
	  hPMCSAngVsMom_2_3_Deg->Fill(p1,0.5*(p2+pU)*(aMCS-aMC)*TMath::DegToRad());
	}
	else if (aAbs >= 3. && aAbs < 10.) {
	  hResMomVertexVsMom_3_10_Deg->Fill(p1,p2-p1);
	  hPDCAVsMom_3_10_Deg->Fill(p1,0.5*(p2+pU)*dca);
	  hPMCSAngVsMom_3_10_Deg->Fill(p1,0.5*(p2+pU)*(aMCS-aMC)*TMath::DegToRad());
	  aMCSMoy += 0.5*(p2+pU)*(aMCS-aMC)*TMath::DegToRad();
	  aMCS2Moy += (0.5*(p2+pU)*(aMCS-aMC)*TMath::DegToRad()) * (0.5*(p2+pU)*(aMCS-aMC)*TMath::DegToRad());
	  dMCSMoy += 0.5*(p2+pU)*(dAbs-pT1/pZ1*AliMUONConstants::AbsZEnd());
	  dMCS2Moy += (0.5*(p2+pU)*(dAbs-pT1/pZ1*AliMUONConstants::AbsZEnd())) * (0.5*(p2+pU)*(dAbs-pT1/pZ1*AliMUONConstants::AbsZEnd()));
	  adMCSMoy += (0.5*(p2+pU)*(aMCS-aMC)*TMath::DegToRad()) * (0.5*(p2+pU)*(dAbs-pT1/pZ1*AliMUONConstants::AbsZEnd()));
	  nMCS++;
	}
	if (aMC < 2.) {
	  hResMomVertexVsMom_0_2_DegMC->Fill(p1,p2-p1);
	  hResMomVertexVsPosAbsEnd_0_2_DegMC->Fill(dAbs,p2-p1);
	  hResSlopeXVertexVsPosAbsEnd_0_2_DegMC->Fill(dAbs,slopex2-slopex1);
	  hResSlopeYVertexVsPosAbsEnd_0_2_DegMC->Fill(dAbs,slopey2-slopey1);
	  hPDCAVsPosAbsEnd_0_2_DegMC->Fill(dAbs,0.5*(p2+pU)*dca);
	  hResEtaVertexVsPosAbsEnd_0_2_DegMC->Fill(dAbs,eta2-eta1);
	  hResPhiVertexVsPosAbsEnd_0_2_DegMC->Fill(dAbs,dPhi);
	}
	else if (aMC >= 2. && aMC < 3) {
	  hResMomVertexVsPosAbsEnd_2_3_DegMC->Fill(dAbs,p2-p1);
	  hResSlopeXVertexVsPosAbsEnd_2_3_DegMC->Fill(dAbs,slopex2-slopex1);
	  hResSlopeYVertexVsPosAbsEnd_2_3_DegMC->Fill(dAbs,slopey2-slopey1);
	  hPDCAVsPosAbsEnd_2_3_DegMC->Fill(dAbs,0.5*(p2+pU)*dca);
	  hResEtaVertexVsPosAbsEnd_2_3_DegMC->Fill(dAbs,eta2-eta1);
	  hResPhiVertexVsPosAbsEnd_2_3_DegMC->Fill(dAbs,dPhi);
	}
	else if (aMC >= 3. && aMC < 10.) {
	  hResMomVertexVsPosAbsEnd_3_10_DegMC->Fill(dAbs,p2-p1);
	  hResSlopeXVertexVsPosAbsEnd_3_10_DegMC->Fill(dAbs,slopex2-slopex1);
	  hResSlopeYVertexVsPosAbsEnd_3_10_DegMC->Fill(dAbs,slopey2-slopey1);
	  hPDCAVsPosAbsEnd_3_10_DegMC->Fill(dAbs,0.5*(p2+pU)*dca);
	  hResEtaVertexVsPosAbsEnd_3_10_DegMC->Fill(dAbs,eta2-eta1);
	  hResPhiVertexVsPosAbsEnd_3_10_DegMC->Fill(dAbs,dPhi);
	}
	hResMomVertexVsAngle->Fill(aAbs,p2-p1);
	hResSlopeXVertexVsAngle->Fill(aAbs,slopex2-slopex1);
	hResSlopeYVertexVsAngle->Fill(aAbs,slopey2-slopey1);
	hPDCAVsAngle->Fill(aAbs,0.5*(p2+pU)*dca);
	hResEtaVertexVsAngle->Fill(aAbs,eta2-eta1);
	hResPhiVertexVsAngle->Fill(aAbs,dPhi);
	hResMomVertexVsMCAngle->Fill(aMC,p2-p1);
	hResSlopeXVertexVsMCAngle->Fill(aMC,slopex2-slopex1);
	hResSlopeYVertexVsMCAngle->Fill(aMC,slopey2-slopey1);
	hPDCAVsMCAngle->Fill(aMC,0.5*(p2+pU)*dca);
	hResEtaVertexVsMCAngle->Fill(aMC,eta2-eta1);
	hResPhiVertexVsMCAngle->Fill(aMC,dPhi);
	
        trackParam = (AliMUONTrackParam*) trackRef->GetTrackParamAtCluster()->First();
        x1 = trackParam->GetNonBendingCoor();
        y1 = trackParam->GetBendingCoor();
        z1 = trackParam->GetZ();
        slopex1 = trackParam->GetNonBendingSlope();
        slopey1 = trackParam->GetBendingSlope();
        pX1 = trackParam->Px();
        pY1 = trackParam->Py();
        pZ1 = trackParam->Pz();
        p1  = trackParam->P();
        pT1 = TMath::Sqrt(pX1*pX1 + pY1*pY1);
	
        trackParam = (AliMUONTrackParam*) trackMatched->GetTrackParamAtCluster()->First();
        x2 = trackParam->GetNonBendingCoor();
        y2 = trackParam->GetBendingCoor();
        z2 = trackParam->GetZ();
        slopex2 = trackParam->GetNonBendingSlope();
        slopey2 = trackParam->GetBendingSlope();
        pX2 = trackParam->Px();
        pY2 = trackParam->Py();
        pZ2 = trackParam->Pz();
        p2  = trackParam->P();
        pT2 = TMath::Sqrt(pX2*pX2 + pY2*pY2);
        
        hResMomFirstCluster->Fill(p2-p1);
	hResMomFirstClusterVsMom->Fill(p1,p2-p1);
	hResPtFirstClusterVsPt->Fill(pT1,pT2-pT1);
	
	hResSlopeXFirstCluster->Fill(slopex2-slopex1);
	hResSlopeYFirstCluster->Fill(slopey2-slopey1);
	hResSlopeXFirstClusterVsMom->Fill(p1,slopex2-slopex1);
	hResSlopeYFirstClusterVsMom->Fill(p1,slopey2-slopey1);
	
	// Fill residuals
	AliMUONTrackParam* trackParamAtCluster1 = (AliMUONTrackParam*) trackMatched->GetTrackParamAtCluster()->First();
	while (trackParamAtCluster1) {
	  AliMUONVCluster* cluster1 = trackParamAtCluster1->GetClusterPtr();
	  AliMUONTrackParam* trackParamAtCluster2 = (AliMUONTrackParam*) trackRef->GetTrackParamAtCluster()->First();
	  while (trackParamAtCluster2) {
	    AliMUONVCluster* cluster2 = trackParamAtCluster2->GetClusterPtr();
	    if (cluster1->GetDetElemId() == cluster2->GetDetElemId()) {
	      hResidualXInCh[cluster1->GetChamberId()]->Fill(cluster1->GetX() - cluster2->GetX());
	      hResidualYInCh[cluster1->GetChamberId()]->Fill(cluster1->GetY() - cluster2->GetY());
	      hResidualXPerDE->Fill(deIndices[cluster1->GetDetElemId()], cluster1->GetX() - cluster2->GetX());
	      hResidualYPerDE->Fill(deIndices[cluster1->GetDetElemId()], cluster1->GetY() - cluster2->GetY());
	      break;
	    }
	    trackParamAtCluster2 = (AliMUONTrackParam*) trackRef->GetTrackParamAtCluster()->After(trackParamAtCluster2);
	  }
	  trackParamAtCluster1 = (AliMUONTrackParam*) trackMatched->GetTrackParamAtCluster()->After(trackParamAtCluster1);
	}
	
      }
      
    } // end loop track ref.

  } // end loop on event  
  cout<<"\rEvent processing... "<<nevents<<" done"<<endl;
  
  // ###################################### compute stuff ###################################### //
  cout<<"\nWhen not specified, resolution at vertex is computed for ";
  if (absorberRegion == 1) cout<<"tracks in the absorber region [2,3] deg."<<endl;
  else if (absorberRegion == 2) cout<<"tracks in the absorber region [3,10] deg."<<endl;
  else cout<<"all tracks"<<endl;
  
  // compute momentum resolution at vertex versus p
  TF1 *f2 = new TF1("f2",langaufun,deltaPAtVtxEdges[0],deltaPAtVtxEdges[1],4);
  Int_t rebinFactorX = TMath::Max(hResMomVertexVsMom->GetNbinsX()/pNBins, 1);
  for (Int_t i = rebinFactorX; i <= hResMomVertexVsMom->GetNbinsX(); i+=rebinFactorX) {
    cout<<"\rFitting momentum residuals at vertex... "<<i/rebinFactorX<<"/"<<pNBins<<flush;
    TH1D *tmp = hResMomVertexVsMom->ProjectionY("tmp",i-rebinFactorX+1,i,"e");
    f2->SetParameters(0.2,0.,(Double_t)tmp->GetEntries(),1.);
    tmp->Fit("f2","WWNQ");
    Double_t fwhm = f2->GetParameter(0);
    Double_t sigma = f2->GetParameter(3);
    Double_t sigmaP = TMath::Sqrt(sigma*sigma + fwhm*fwhm/(8.*log(2.)));
    Int_t rebin = TMath::Max(Int_t(0.5*sigmaP/tmp->GetBinWidth(1)),1);
    while (deltaPAtVtxNBins%rebin!=0) rebin--;
    tmp->Rebin(rebin);
    tmp->Fit("f2","NQ");
    fwhm = f2->GetParameter(0);
    sigma = f2->GetParameter(3);
    sigmaP = TMath::Sqrt(sigma*sigma + fwhm*fwhm/(8.*log(2.)));
    Double_t fwhmErr = f2->GetParError(0);
    Double_t sigmaErr = f2->GetParError(3);
    Double_t sigmaPErr = TMath::Sqrt(sigma*sigma*sigmaErr*sigmaErr + fwhm*fwhm*fwhmErr*fwhmErr/(64.*log(2.)*log(2.))) / sigmaP;
    hResMomVertexVsMom->GetXaxis()->SetRange(i-rebinFactorX+1,i);
    Double_t p = (tmp->GetEntries() > 0) ? hResMomVertexVsMom->GetMean() : 0.5 * (hResMomVertexVsMom->GetBinLowEdge(i-rebinFactorX+1) + hResMomVertexVsMom->GetBinLowEdge(i+1));
    hResMomVertexVsMom->GetXaxis()->SetRange();
    Double_t pErr[2] = {p-hResMomVertexVsMom->GetBinLowEdge(i-rebinFactorX+1), hResMomVertexVsMom->GetBinLowEdge(i+1)-p};
    gMeanResMomVertexVsMom->SetPoint(i/rebinFactorX-1, p, tmp->GetMean());
    gMeanResMomVertexVsMom->SetPointError(i/rebinFactorX-1, pErr[0], pErr[1], tmp->GetMeanError(), tmp->GetMeanError());
    gMostProbResMomVertexVsMom->SetPoint(i/rebinFactorX-1, p, -f2->GetParameter(1));
    gMostProbResMomVertexVsMom->SetPointError(i/rebinFactorX-1, pErr[0], pErr[1], f2->GetParError(1), f2->GetParError(1));
    gSigmaResMomVertexVsMom->SetPoint(i/rebinFactorX-1, p, 100.*sigmaP/p);
    gSigmaResMomVertexVsMom->SetPointError(i/rebinFactorX-1, pErr[0], pErr[1], 100.*sigmaPErr/p, 100.*sigmaPErr/p);
    delete tmp;
  }
  cout<<"\rFitting momentum residuals at vertex... "<<pNBins<<"/"<<pNBins<<endl;
  
  // compute momentum relative resolution at first cluster versus p
  FitGausResVsMom(hResMomFirstClusterVsMom, pNBins, 0., 1., "momentum residuals at first cluster", gMeanResMomFirstClusterVsMom, gSigmaResMomFirstClusterVsMom);
  rebinFactorX = TMath::Max(hResMomFirstClusterVsMom->GetNbinsX()/pNBins, 1);
  for (Int_t i = rebinFactorX; i <= hResMomFirstClusterVsMom->GetNbinsX(); i+=rebinFactorX) {
    Double_t x,y;
    gSigmaResMomFirstClusterVsMom->GetPoint(i/rebinFactorX-1, x, y);
    gSigmaResMomFirstClusterVsMom->SetPoint(i/rebinFactorX-1, x, 100.*y/x);
    gSigmaResMomFirstClusterVsMom->SetPointEYlow(i/rebinFactorX-1, 100.*gSigmaResMomFirstClusterVsMom->GetErrorYlow(i/rebinFactorX-1)/x);
    gSigmaResMomFirstClusterVsMom->SetPointEYhigh(i/rebinFactorX-1, 100.*gSigmaResMomFirstClusterVsMom->GetErrorYhigh(i/rebinFactorX-1)/x);
  }
  
  // compute slopeX resolution at vertex versus p
  FitGausResVsMom(hResSlopeXVertexVsMom, pNBins, 0., 2.e-3, "slopeX residuals at vertex", gMeanResSlopeXVertexVsMom, gSigmaResSlopeXVertexVsMom);
  
  // compute slopeY resolution at vertex versus p
  FitGausResVsMom(hResSlopeYVertexVsMom, pNBins, 0., 2.e-3, "slopeY residuals at vertex", gMeanResSlopeYVertexVsMom, gSigmaResSlopeYVertexVsMom);
  
  // compute slopeX resolution at first cluster versus p
  FitGausResVsMom(hResSlopeXFirstClusterVsMom, pNBins, 0., 3.e-4, "slopeX residuals at first cluster", gMeanResSlopeXFirstClusterVsMom, gSigmaResSlopeXFirstClusterVsMom);
  
  // compute slopeY resolution at first cluster versus p
  FitGausResVsMom(hResSlopeYFirstClusterVsMom, pNBins, 0., 2.e-4, "slopeY residuals at first cluster", gMeanResSlopeYFirstClusterVsMom, gSigmaResSlopeYFirstClusterVsMom);
  
  // compute p*DCA resolution in the region [2,3] deg at absorber end
  FitPDCAVsMom(hPDCAVsMom_2_3_Deg, pNBins, "p*DCA (tracks in [2,3] deg.)", gMeanPDCAVsMom_2_3_Deg, gSigmaPDCAVsMom_2_3_Deg);
  
  // compute p*DCA resolution in the region [3,10] deg at absorber end
  FitPDCAVsMom(hPDCAVsMom_3_10_Deg, pNBins, "p*DCA (tracks in [3,10] deg.)", gMeanPDCAVsMom_3_10_Deg, gSigmaPDCAVsMom_3_10_Deg);
  
  // compute MCS angular dispersion in the region [2,3] deg at absorber end
  FitGausResVsMom(hPMCSAngVsMom_2_3_Deg, pNBins, 0., 2.e-3, "p*MCSAngle (tracks in [2,3] deg.)", gMeanPMCSAngVsMom_2_3_Deg, gSigmaPMCSAngVsMom_2_3_Deg);
  
  // compute MCS angular dispersion in the region [3,10] deg at absorber end
  FitGausResVsMom(hPMCSAngVsMom_3_10_Deg, pNBins, 0., 2.e-3, "p*MCSAngle (tracks in [3,10] deg.)", gMeanPMCSAngVsMom_3_10_Deg, gSigmaPMCSAngVsMom_3_10_Deg);
  
  // compute eta resolution at vertex versus p
  FitGausResVsMom(hResEtaVertexVsMom, pNBins, 0., 0.1, "eta residuals at vertex", gMeanResEtaVertexVsMom, gSigmaResEtaVertexVsMom);
  
  // compute phi resolution at vertex versus p
  FitGausResVsMom(hResPhiVertexVsMom, pNBins, 0., 0.01, "phi residuals at vertex", gMeanResPhiVertexVsMom, gSigmaResPhiVertexVsMom);
  
  // compute cluster-track residual mean and dispersion per chamber
  Double_t clusterResPerCh[10][2];
  for (Int_t i = 0; i < AliMUONConstants::NTrackingCh(); i++) {
    FillResidual(hResidualXInCh[i], i, clusterResPerCh[i][0], gResidualXPerChMean, gResidualXPerChSigma, kTRUE, fitResiduals);
    FillResidual(hResidualYInCh[i], i, clusterResPerCh[i][1], gResidualYPerChMean, gResidualYPerChSigma, kTRUE, fitResiduals);
  }
  
  // compute cluster-track residual mean and dispersion per DE
  Double_t clusterResPerDE[200][2];
  for (Int_t i = 0; i < deNBins; i++) {
    TH1D *tmp = hResidualXPerDE->ProjectionY("tmp",i+1,i+1,"e");
    FillResidual(tmp, i, clusterResPerDE[i][0], gResidualXPerDEMean, gResidualXPerDESigma, kTRUE, fitResiduals);
    delete tmp;
    tmp = hResidualYPerDE->ProjectionY("tmp",i+1,i+1,"e");
    FillResidual(tmp, i, clusterResPerDE[i][1], gResidualYPerDEMean, gResidualYPerDESigma, kTRUE, fitResiduals);
    delete tmp;
  }
  
  // ###################################### display histograms ###################################### //
  // diplay momentum residuals
  TCanvas* cResMom = DrawVsAng("cResMom", "momentum residual at vertex in 3 angular regions", hResMomVertex, hResMomVertexVsAngle);
  TCanvas* cResMomMC = DrawVsAng("cResMomMC", "momentum residual at vertex in 3 MC angular regions", hResMomVertex, hResMomVertexVsMCAngle);
  TCanvas* cResMomVsPos = DrawVsPos("cResMomVsPos", "momentum residual at vertex versus position at absorber end in 3 MC angular regions",
				    hResMomVertexVsPosAbsEnd_0_2_DegMC, hResMomVertexVsPosAbsEnd_2_3_DegMC, hResMomVertexVsPosAbsEnd_3_10_DegMC);
  TCanvas* cResMom_2_3_Deg = DrawResMomVsMom("cResMom_2_3_Deg", "momentum residual for tracks between 2 and 3 degrees",
					     hResMomVertexVsMom_2_3_Deg, 10, f2, "momentum residuals at vertex (tracks in [2,3] deg.)");
  TCanvas* cResMom_3_10_Deg = DrawResMomVsMom("cResMom_3_10_Deg", "momentum residual for tracks between 3 and 10 degrees",
					      hResMomVertexVsMom_3_10_Deg, 10, f2, "momentum residuals at vertex (tracks in [3,10] deg.)");
  TCanvas* cResMom_0_2_DegMC = DrawResMomVsMom("cResMom_0_2_DegMC", "momentum residuals for tracks with MC angle < 2 degrees", hResMomVertexVsMom_0_2_DegMC, 5);
  
  // diplay slopeX residuals
  TCanvas* cResSlopeX = DrawVsAng("cResSlopeX", "slope_{X} residual at vertex in 3 angular regions", hResSlopeXVertex, hResSlopeXVertexVsAngle);
  TCanvas* cResSlopeXMC = DrawVsAng("cResSlopeXMC", "slope_{X} residual at vertex in 3 MC angular regions", hResSlopeXVertex, hResSlopeXVertexVsMCAngle);
  TCanvas* cResSlopeXVsPos = DrawVsPos("cResSlopeXVsPos", "slope_{X} residual at vertex versus position at absorber end in 3 MC angular regions",
				       hResSlopeXVertexVsPosAbsEnd_0_2_DegMC, hResSlopeXVertexVsPosAbsEnd_2_3_DegMC, hResSlopeXVertexVsPosAbsEnd_3_10_DegMC);
  
  // diplay slopeY residuals
  TCanvas* cResSlopeY = DrawVsAng("cResSlopeY", "slope_{Y} residual at vertex in 3 angular regions", hResSlopeYVertex, hResSlopeYVertexVsAngle);
  TCanvas* cResSlopeYMC = DrawVsAng("cResSlopeYMC", "slope_{Y} residual at vertex in 3 MC angular regions", hResSlopeYVertex, hResSlopeYVertexVsMCAngle);
  TCanvas* cResSlopeYVsPos = DrawVsPos("cResSlopeYVsPos", "slope_{Y} residual at vertex versus position at absorber end in 3 MC angular regions",
				       hResSlopeYVertexVsPosAbsEnd_0_2_DegMC, hResSlopeYVertexVsPosAbsEnd_2_3_DegMC, hResSlopeYVertexVsPosAbsEnd_3_10_DegMC);
  
  // diplay P*DCA
  TCanvas* cPDCA = DrawVsAng("cPDCA", "p #times DCA in 3 angular regions", hPDCA, hPDCAVsAngle);
  TCanvas* cPDCAMC = DrawVsAng("cPDCAMC", "p #times DCA in 3 MC angular regions", hPDCA, hPDCAVsMCAngle);
  TCanvas* cPDCAVsPos = DrawVsPos("cPDCAVsPos", "p #times DCA versus position at absorber end in 3 MC angular regions",
				  hPDCAVsPosAbsEnd_0_2_DegMC, hPDCAVsPosAbsEnd_2_3_DegMC, hPDCAVsPosAbsEnd_3_10_DegMC);
  
  // diplay eta residuals
  TCanvas* cResEta = DrawVsAng("cResEta", "eta residual at vertex in 3 angular regions", hResEtaVertex, hResEtaVertexVsAngle);
  TCanvas* cResEtaMC = DrawVsAng("cResEtaMC", "eta residual at vertex in 3 MC angular regions", hResEtaVertex, hResEtaVertexVsMCAngle);
  TCanvas* cResEtaVsPos = DrawVsPos("cResEtaVsPos", "eta residual at vertex versus position at absorber end in 3 MC angular regions",
				    hResEtaVertexVsPosAbsEnd_0_2_DegMC, hResEtaVertexVsPosAbsEnd_2_3_DegMC, hResEtaVertexVsPosAbsEnd_3_10_DegMC);
  
  // diplay phi residuals
  TCanvas* cResPhi = DrawVsAng("cResPhi", "phi residual at vertex in 3 angular regions", hResPhiVertex, hResPhiVertexVsAngle);
  TCanvas* cResPhiMC = DrawVsAng("cResPhiMC", "phi residual at vertex in 3 MC angular regions", hResPhiVertex, hResPhiVertexVsMCAngle);
  TCanvas* cResPhiVsPos = DrawVsPos("cResPhiVsPos", "phi residual at vertex versus position at absorber end in 3 MC angular regions",
				    hResPhiVertexVsPosAbsEnd_0_2_DegMC, hResPhiVertexVsPosAbsEnd_2_3_DegMC, hResPhiVertexVsPosAbsEnd_3_10_DegMC);
  
  // ###################################### save histogram ###################################### //
  histoFile->Write();
  
  histoFile->cd("momentumAtVertex");
  gMeanResMomVertexVsMom->Write();
  gMostProbResMomVertexVsMom->Write();
  gSigmaResMomVertexVsMom->Write();
  cResMom->Write();
  cResMomMC->Write();
  cResMomVsPos->Write();
  cResMom_2_3_Deg->Write();
  cResMom_3_10_Deg->Write();
  cResMom_0_2_DegMC->Write();
  
  histoFile->cd("slopesAtVertex");
  gMeanResSlopeXVertexVsMom->Write();
  gMeanResSlopeYVertexVsMom->Write();
  gSigmaResSlopeXVertexVsMom->Write();
  gSigmaResSlopeYVertexVsMom->Write();
  cResSlopeX->Write();
  cResSlopeY->Write();
  cResSlopeXMC->Write();
  cResSlopeYMC->Write();
  cResSlopeXVsPos->Write();
  cResSlopeYVsPos->Write();
  
  histoFile->cd("DCA");
  gMeanPDCAVsMom_2_3_Deg->Write();
  gSigmaPDCAVsMom_2_3_Deg->Write();
  gMeanPDCAVsMom_3_10_Deg->Write();
  gSigmaPDCAVsMom_3_10_Deg->Write();
  gMeanPMCSAngVsMom_2_3_Deg->Write();
  gSigmaPMCSAngVsMom_2_3_Deg->Write();
  gMeanPMCSAngVsMom_3_10_Deg->Write();
  gSigmaPMCSAngVsMom_3_10_Deg->Write();
  cPDCA->Write();
  cPDCAMC->Write();
  cPDCAVsPos->Write();
  
  histoFile->cd("etaAtVertex");
  gMeanResEtaVertexVsMom->Write();
  gSigmaResEtaVertexVsMom->Write();
  cResEta->Write();
  cResEtaMC->Write();
  cResEtaVsPos->Write();
  
  histoFile->cd("phiAtVertex");
  gMeanResPhiVertexVsMom->Write();
  gSigmaResPhiVertexVsMom->Write();
  cResPhi->Write();
  cResPhiMC->Write();
  cResPhiVsPos->Write();
  
  histoFile->cd("momentumAtFirstCluster");
  gMeanResMomFirstClusterVsMom->Write();
  gSigmaResMomFirstClusterVsMom->Write();
  
  histoFile->cd("slopesAtFirstCluster");
  gMeanResSlopeXFirstClusterVsMom->Write();
  gMeanResSlopeYFirstClusterVsMom->Write();
  gSigmaResSlopeXFirstClusterVsMom->Write();
  gSigmaResSlopeYFirstClusterVsMom->Write();
  
  histoFile->cd("clusters");
  gResidualXPerChMean->Write();
  gResidualXPerChSigma->Write();
  gResidualYPerChMean->Write();
  gResidualYPerChSigma->Write();
  gResidualXPerDEMean->Write();
  gResidualXPerDESigma->Write();
  gResidualYPerDEMean->Write();
  gResidualYPerDESigma->Write();
  
  histoFile->Close();
  
  // ###################################### clean memory ###################################### //
  delete cResMom;
  delete cResMomMC;
  delete cResMomVsPos;
  delete cResMom_2_3_Deg;
  delete cResMom_3_10_Deg;
  delete cResMom_0_2_DegMC;
  delete cResSlopeX;
  delete cResSlopeY;
  delete cResSlopeXMC;
  delete cResSlopeYMC;
  delete cResSlopeXVsPos;
  delete cResSlopeYVsPos;
  delete cPDCA;
  delete cPDCAMC;
  delete cPDCAVsPos;
  delete cResEta;
  delete cResEtaMC;
  delete cResEtaVsPos;
  delete cResPhi;
  delete cResPhiMC;
  delete cResPhiVsPos;
  
  // ###################################### print statistics ###################################### //
  printf("\n");
  printf("nb of reconstructible tracks: %d \n", nReconstructibleTracks);
  printf("nb of reconstructed tracks: %d \n", nReconstructedTracks);
  printf("nb of reconstructible tracks which are reconstructed: %d \n", nReconstructibleTracksCheck);
  
  aMCSMoy /= (Double_t) nMCS;
  aMCS2Moy /= (Double_t) nMCS;
  dMCSMoy /= (Double_t) nMCS;
  dMCS2Moy /= (Double_t) nMCS;
  adMCSMoy /= (Double_t) nMCS;
  Double_t sigma2_ThetaMCS = aMCS2Moy - aMCSMoy*aMCSMoy;
  Double_t sigma2_PosMCS = dMCS2Moy - dMCSMoy*dMCSMoy;
  Double_t cov_ThetaPosMCS = - (adMCSMoy - aMCSMoy*dMCSMoy);
  printf("\nmultiple scattering of tracks between 3 and 10 deg. at absorber end:\n");
  printf(" sigma_ThetaMCS = %f\n", TMath::Sqrt(sigma2_ThetaMCS));
  printf(" sigma_PosMCS = %f\n", TMath::Sqrt(sigma2_PosMCS));
  printf(" cov_ThetaPosMCS = %f\n", cov_ThetaPosMCS);
  printf(" --> sigma_DCA = %f\n", TMath::Sqrt(AliMUONConstants::AbsZEnd()*AliMUONConstants::AbsZEnd()*sigma2_ThetaMCS
					      - 2.*AliMUONConstants::AbsZEnd()*cov_ThetaPosMCS + sigma2_PosMCS));
  
  printf("\nchamber resolution:\n");
  printf(" - non-bending:");
  for (Int_t i = 0; i < AliMUONConstants::NTrackingCh(); i++) printf((i==0)?" %5.3f":", %5.3f",clusterResPerCh[i][0]);
  printf("\n -     bending:");
  for (Int_t i = 0; i < AliMUONConstants::NTrackingCh(); i++) printf((i==0)?" %6.4f":", %6.4f",clusterResPerCh[i][1]);
  printf("\n\n");
  
  printf("\nDE resolution:\n");
  printf(" - non-bending:");
  for (Int_t i = 0; i < deNBins; i++) printf((i==0)?" %5.3f":", %5.3f",clusterResPerDE[i][0]);
  printf("\n -     bending:");
  for (Int_t i = 0; i < deNBins; i++) printf((i==0)?" %6.4f":", %6.4f",clusterResPerDE[i][1]);
  printf("\n\n");
}

//------------------------------------------------------------------------------------
Double_t langaufun(Double_t *x, Double_t *par) {
  
  //Fit parameters:
  //par[0]=Width (scale) parameter of Landau density
  //par[1]=Most Probable (MP, location) parameter of Landau density
  //par[2]=Total area (integral -inf to inf, normalization constant)
  //par[3]=Width (sigma) of convoluted Gaussian function
  //
  //In the Landau distribution (represented by the CERNLIB approximation), 
  //the maximum is located at x=-0.22278298 with the location parameter=0.
  //This shift is corrected within this function, so that the actual
  //maximum is identical to the MP parameter.
  
  // Numeric constants
  Double_t invsq2pi = 0.3989422804014;   // (2 pi)^(-1/2)
  Double_t mpshift  = -0.22278298;       // Landau maximum location
  
  // Control constants
  Double_t np = 100.0; // number of convolution steps
  Double_t sc = 5.0;   // convolution extends to +-sc Gaussian sigmas
  
  // Variables
  Double_t xx;
  Double_t mpc;
  Double_t fland;
  Double_t sum = 0.0;
  Double_t xlow,xupp;
  Double_t step;
  Double_t i;
  
  
  // MP shift correction
  mpc = par[1] - mpshift * par[0]; 
  
  // Range of convolution integral
  xlow = x[0] - sc * par[3];
  xupp = x[0] + sc * par[3];
  
  step = (xupp-xlow) / np;
  
  // Convolution integral of Landau and Gaussian by sum
  for(i=1.0; i<=np/2; i++) {
    xx = xlow + (i-.5) * step;
    //change x -> -x because the tail of the Landau is at the left here...
    fland = TMath::Landau(-xx,mpc,par[0]) / par[0];
    sum += fland * TMath::Gaus(x[0],xx,par[3]);
    
    //change x -> -x because the tail of the Landau is at the left here...
    xx = xupp - (i-.5) * step;
    fland = TMath::Landau(-xx,mpc,par[0]) / par[0];
    sum += fland * TMath::Gaus(x[0],xx,par[3]);
  }
  
  return (par[2] * step * sum * invsq2pi / par[3]);
}

//------------------------------------------------------------------------------------
void FitGausResVsMom(TH2* h, Int_t nBins, const Double_t mean0, const Double_t sigma0,
		     const char* fitting, TGraphAsymmErrors* gMean, TGraphAsymmErrors* gSigma)
{
  /// generic function to fit residuals versus momentum with a gaussian
  static TF1* fGaus = 0x0;
  if (!fGaus) fGaus = new TF1("fGaus","gaus");
  
  Int_t rebinFactorX = TMath::Max(h->GetNbinsX()/nBins, 1);
  for (Int_t i = rebinFactorX; i <= h->GetNbinsX(); i+=rebinFactorX) {
    cout<<Form("\rFitting %s... %d/%d",fitting,i/rebinFactorX,nBins)<<flush;
    TH1D *tmp = h->ProjectionY("tmp",i-rebinFactorX+1,i,"e");
    fGaus->SetParameters(tmp->GetEntries(), mean0, sigma0);
    tmp->Fit("fGaus","WWNQ");
    Int_t rebin = TMath::Max(Int_t(0.5*fGaus->GetParameter(2)/tmp->GetBinWidth(1)),1);
    while (tmp->GetNbinsX()%rebin!=0) rebin--;
    tmp->Rebin(rebin);
    tmp->Fit("fGaus","NQ");
    h->GetXaxis()->SetRange(i-rebinFactorX+1,i);
    Double_t p = (tmp->GetEntries() > 0) ? h->GetMean() : 0.5 * (h->GetBinLowEdge(i-rebinFactorX+1) + h->GetBinLowEdge(i+1));
    h->GetXaxis()->SetRange();
    Double_t pErr[2] = {p-h->GetBinLowEdge(i-rebinFactorX+1), h->GetBinLowEdge(i+1)-p};
    gMean->SetPoint(i/rebinFactorX-1, p, fGaus->GetParameter(1));
    gMean->SetPointError(i/rebinFactorX-1, pErr[0], pErr[1], fGaus->GetParError(1), fGaus->GetParError(1));
    gSigma->SetPoint(i/rebinFactorX-1, p, fGaus->GetParameter(2));
    gSigma->SetPointError(i/rebinFactorX-1, pErr[0], pErr[1], fGaus->GetParError(2), fGaus->GetParError(2));
    delete tmp;
  }
  cout<<Form("\rFitting %s... %d/%d",fitting,nBins,nBins)<<endl;
}

//------------------------------------------------------------------------------------
void FitPDCAVsMom(TH2* h, Int_t nBins, const char* fitting, TGraphAsymmErrors* gMean, TGraphAsymmErrors* gSigma)
{
  /// generic function to fit p*DCA distributions
  static TF1* fPGaus = 0x0;
  if (!fPGaus) fPGaus = new TF1("fPGaus","x*gaus");
  
  Int_t rebinFactorX = TMath::Max(h->GetNbinsX()/nBins, 1);
  for (Int_t i = rebinFactorX; i <= h->GetNbinsX(); i+=rebinFactorX) {
    cout<<Form("\rFitting %s... %d/%d",fitting,i/rebinFactorX,nBins)<<flush;
    TH1D *tmp = h->ProjectionY("tmp",i-rebinFactorX+1,i,"e");
    fPGaus->SetParameters(1.,0.,80.);
    Int_t rebin = 25.*(tmp->GetNbinsX()/(tmp->GetBinLowEdge(tmp->GetNbinsX()+1)-tmp->GetBinLowEdge(1)));
    while (tmp->GetNbinsX()%rebin!=0) rebin--;
    tmp->Rebin(rebin);
    tmp->Fit("fPGaus","NQ");
    h->GetXaxis()->SetRange(i-rebinFactorX+1,i);
    Double_t p = (tmp->GetEntries() > 0) ? h->GetMean() : 0.5 * (h->GetBinLowEdge(i-rebinFactorX+1) + h->GetBinLowEdge(i+1));
    h->GetXaxis()->SetRange();
    Double_t pErr[2] = {p-h->GetBinLowEdge(i-rebinFactorX+1), h->GetBinLowEdge(i+1)-p};
    gMean->SetPoint(i/rebinFactorX-1, p, fPGaus->GetParameter(1));
    gMean->SetPointError(i/rebinFactorX-1, pErr[0], pErr[1], fPGaus->GetParError(1), fPGaus->GetParError(1));
    gSigma->SetPoint(i/rebinFactorX-1, p, fPGaus->GetParameter(2));
    gSigma->SetPointError(i/rebinFactorX-1, pErr[0], pErr[1], fPGaus->GetParError(2), fPGaus->GetParError(2));
    delete tmp;
  }
  cout<<Form("\rFitting %s... %d/%d",fitting,nBins,nBins)<<endl;
}

//------------------------------------------------------------------------------------
void FillResidual(TH1* h, Int_t i, Double_t& sigma, TGraphErrors* gMean, TGraphErrors* gSigma, Bool_t correctForSystematics, Bool_t fitResiduals)
{
  /// fill graphs with residual mean and sigma
  static TF1* fRGaus = 0x0;
  Double_t mean, meanErr, sigmaErr;
  
  if (fitResiduals) {
    
    if (!fRGaus) fRGaus = new TF1("fRGaus","gaus");
    fRGaus->SetParameters(h->GetEntries(), 0., 0.1);
    h->Fit("fRGaus", "WWNQ");
    mean = fRGaus->GetParameter(1);
    meanErr = fRGaus->GetParError(1);
    sigma = fRGaus->GetParameter(2);
    sigmaErr = fRGaus->GetParError(2);
    
  } else {
    
    Zoom(h);
    mean = h->GetMean();
    meanErr = h->GetMeanError();
    sigma = h->GetRMS();
    sigmaErr = h->GetRMSError();
    h->GetXaxis()->SetRange(0,0);
    
  }
  
  gMean->SetPoint(i, i+1, mean);
  gMean->SetPointError(i, 0., meanErr);
  if (correctForSystematics) {
    Double_t s = TMath::Sqrt(mean*mean + sigma*sigma);
    sigmaErr = (s>0.) ? TMath::Sqrt(sigma*sigma*sigmaErr*sigmaErr + mean*mean*meanErr*meanErr) / s : 0.;
    sigma = s;
  }
  gSigma->SetPoint(i, i+1, sigma);
  gSigma->SetPointError(i, 0., sigmaErr);
}

//------------------------------------------------------------------------------------
TCanvas* DrawVsAng(const char* name, const char* title, TH1* h1, TH2* h2)
{
  /// generic function to draw histograms versus absorber angular region
  TCanvas* c = new TCanvas(name, title);
  c->cd();
  h1->Draw();
  TH1D *proj1 = h2->ProjectionY(Form("%s_proj_0_2",h2->GetName()),1,2);
  proj1->Draw("sames");
  proj1->SetLineColor(2);
  TH1D *proj2 = h2->ProjectionY(Form("%s_proj_2_3",h2->GetName()),3,3);
  proj2->Draw("sames");
  proj2->SetLineColor(4);
  TH1D *proj3 = h2->ProjectionY(Form("%s__proj_3_10",h2->GetName()),4,10);
  proj3->Draw("sames");
  proj3->SetLineColor(3);
  return c;
}

//------------------------------------------------------------------------------------
TCanvas* DrawVsPos(const char* name, const char* title, TH2* h1, TH2* h2, TH2* h3)
{
  /// generic function to draw histograms versus position at absorber end
  TCanvas* c = new TCanvas(name, title);
  c->cd();
  h1->Draw();
  h1->SetMarkerColor(2);
  h2->Draw("sames");
  h2->SetMarkerColor(4);
  h3->Draw("sames");
  h3->SetMarkerColor(3);
  return c;
}

//------------------------------------------------------------------------------------
TCanvas* DrawResMomVsMom(const char* name, const char* title, TH2* h, Int_t nBins, TF1* f2, const char* fitting)
{
  /// generic function to draw and eventually fit momentum residuals versus momentum
  TLegend* l = new TLegend(0.15,0.25,0.3,0.85);
  TCanvas* c = new TCanvas(name, title);
  c->cd();
  TH1D* proj = 0x0;
  h->Sumw2();
  Int_t rebinFactorX = TMath::Max(h->GetNbinsX()/nBins, 1);
  for (Int_t i = rebinFactorX; i <= h->GetNbinsX(); i+=rebinFactorX) {
    if (f2) cout<<Form("\rFitting %s... %d/%d",fitting,i/rebinFactorX,nBins)<<flush;
    proj = h->ProjectionY(Form("%s_%d",h->GetName(),i/rebinFactorX),i-rebinFactorX+1,i);
    if (proj->GetEntries() > 0) proj->Scale(1./proj->GetEntries());
    proj->Draw((i==rebinFactorX)?"hist":"histsames");
    proj->SetLineColor(i/rebinFactorX);
    if (f2) {
      f2->SetParameters(0.2,0.,1.,1.);
      f2->SetLineColor(i/rebinFactorX);
      proj->Fit("f2","WWNQ","sames");
      Double_t fwhm = f2->GetParameter(0);
      Double_t sigma = f2->GetParameter(3);
      Double_t sigmaP = TMath::Sqrt(sigma*sigma + fwhm*fwhm/(8.*log(2.)));
      Int_t rebin = TMath::Max(Int_t(0.5*sigmaP/proj->GetBinWidth(1)),1);
      while (proj->GetNbinsX()%rebin!=0) rebin--;
      proj->Rebin(rebin);
      proj->Scale(1./rebin);
      proj->Fit("f2","Q","sames");
    } else proj->SetLineWidth(2);
    Double_t p = 0.5 * (h->GetBinLowEdge(i-rebinFactorX+1) + h->GetBinLowEdge(i+1));
    l->AddEntry(proj,Form("%5.1f GeV",p));
  }
  if (f2) cout<<Form("\rFitting %s... %d/%d",fitting,nBins,nBins)<<endl;
  l->Draw("same");
  return c;
}

//------------------------------------------------------------------------------------
void Zoom(TH1* h, Double_t fractionCut)
{
  /// Reduce the range of the histogram by removing a given fration of the statistic at each edge
  Int_t maxEventsCut = fractionCut * h->GetEntries();
  Int_t nBins = h->GetNbinsX();
  
  // set low edge  
  Int_t minBin;
  Int_t eventsCut = 0;
  for (minBin = 1; minBin <= nBins; minBin++) {
    eventsCut += h->GetBinContent(minBin);
    if (eventsCut > maxEventsCut) break;
  }
  
  // set high edge
  Int_t maxBin;
  eventsCut = 0;
  for (maxBin = nBins; maxBin >= 1; maxBin--) {
    eventsCut += h->GetBinContent(maxBin);
    if (eventsCut > maxEventsCut) break;
  }
  
  // set new axis range
  h->GetXaxis()->SetRange(--minBin, ++maxBin);
}

 MUONRecoCheck.C:1
 MUONRecoCheck.C:2
 MUONRecoCheck.C:3
 MUONRecoCheck.C:4
 MUONRecoCheck.C:5
 MUONRecoCheck.C:6
 MUONRecoCheck.C:7
 MUONRecoCheck.C:8
 MUONRecoCheck.C:9
 MUONRecoCheck.C:10
 MUONRecoCheck.C:11
 MUONRecoCheck.C:12
 MUONRecoCheck.C:13
 MUONRecoCheck.C:14
 MUONRecoCheck.C:15
 MUONRecoCheck.C:16
 MUONRecoCheck.C:17
 MUONRecoCheck.C:18
 MUONRecoCheck.C:19
 MUONRecoCheck.C:20
 MUONRecoCheck.C:21
 MUONRecoCheck.C:22
 MUONRecoCheck.C:23
 MUONRecoCheck.C:24
 MUONRecoCheck.C:25
 MUONRecoCheck.C:26
 MUONRecoCheck.C:27
 MUONRecoCheck.C:28
 MUONRecoCheck.C:29
 MUONRecoCheck.C:30
 MUONRecoCheck.C:31
 MUONRecoCheck.C:32
 MUONRecoCheck.C:33
 MUONRecoCheck.C:34
 MUONRecoCheck.C:35
 MUONRecoCheck.C:36
 MUONRecoCheck.C:37
 MUONRecoCheck.C:38
 MUONRecoCheck.C:39
 MUONRecoCheck.C:40
 MUONRecoCheck.C:41
 MUONRecoCheck.C:42
 MUONRecoCheck.C:43
 MUONRecoCheck.C:44
 MUONRecoCheck.C:45
 MUONRecoCheck.C:46
 MUONRecoCheck.C:47
 MUONRecoCheck.C:48
 MUONRecoCheck.C:49
 MUONRecoCheck.C:50
 MUONRecoCheck.C:51
 MUONRecoCheck.C:52
 MUONRecoCheck.C:53
 MUONRecoCheck.C:54
 MUONRecoCheck.C:55
 MUONRecoCheck.C:56
 MUONRecoCheck.C:57
 MUONRecoCheck.C:58
 MUONRecoCheck.C:59
 MUONRecoCheck.C:60
 MUONRecoCheck.C:61
 MUONRecoCheck.C:62
 MUONRecoCheck.C:63
 MUONRecoCheck.C:64
 MUONRecoCheck.C:65
 MUONRecoCheck.C:66
 MUONRecoCheck.C:67
 MUONRecoCheck.C:68
 MUONRecoCheck.C:69
 MUONRecoCheck.C:70
 MUONRecoCheck.C:71
 MUONRecoCheck.C:72
 MUONRecoCheck.C:73
 MUONRecoCheck.C:74
 MUONRecoCheck.C:75
 MUONRecoCheck.C:76
 MUONRecoCheck.C:77
 MUONRecoCheck.C:78
 MUONRecoCheck.C:79
 MUONRecoCheck.C:80
 MUONRecoCheck.C:81
 MUONRecoCheck.C:82
 MUONRecoCheck.C:83
 MUONRecoCheck.C:84
 MUONRecoCheck.C:85
 MUONRecoCheck.C:86
 MUONRecoCheck.C:87
 MUONRecoCheck.C:88
 MUONRecoCheck.C:89
 MUONRecoCheck.C:90
 MUONRecoCheck.C:91
 MUONRecoCheck.C:92
 MUONRecoCheck.C:93
 MUONRecoCheck.C:94
 MUONRecoCheck.C:95
 MUONRecoCheck.C:96
 MUONRecoCheck.C:97
 MUONRecoCheck.C:98
 MUONRecoCheck.C:99
 MUONRecoCheck.C:100
 MUONRecoCheck.C:101
 MUONRecoCheck.C:102
 MUONRecoCheck.C:103
 MUONRecoCheck.C:104
 MUONRecoCheck.C:105
 MUONRecoCheck.C:106
 MUONRecoCheck.C:107
 MUONRecoCheck.C:108
 MUONRecoCheck.C:109
 MUONRecoCheck.C:110
 MUONRecoCheck.C:111
 MUONRecoCheck.C:112
 MUONRecoCheck.C:113
 MUONRecoCheck.C:114
 MUONRecoCheck.C:115
 MUONRecoCheck.C:116
 MUONRecoCheck.C:117
 MUONRecoCheck.C:118
 MUONRecoCheck.C:119
 MUONRecoCheck.C:120
 MUONRecoCheck.C:121
 MUONRecoCheck.C:122
 MUONRecoCheck.C:123
 MUONRecoCheck.C:124
 MUONRecoCheck.C:125
 MUONRecoCheck.C:126
 MUONRecoCheck.C:127
 MUONRecoCheck.C:128
 MUONRecoCheck.C:129
 MUONRecoCheck.C:130
 MUONRecoCheck.C:131
 MUONRecoCheck.C:132
 MUONRecoCheck.C:133
 MUONRecoCheck.C:134
 MUONRecoCheck.C:135
 MUONRecoCheck.C:136
 MUONRecoCheck.C:137
 MUONRecoCheck.C:138
 MUONRecoCheck.C:139
 MUONRecoCheck.C:140
 MUONRecoCheck.C:141
 MUONRecoCheck.C:142
 MUONRecoCheck.C:143
 MUONRecoCheck.C:144
 MUONRecoCheck.C:145
 MUONRecoCheck.C:146
 MUONRecoCheck.C:147
 MUONRecoCheck.C:148
 MUONRecoCheck.C:149
 MUONRecoCheck.C:150
 MUONRecoCheck.C:151
 MUONRecoCheck.C:152
 MUONRecoCheck.C:153
 MUONRecoCheck.C:154
 MUONRecoCheck.C:155
 MUONRecoCheck.C:156
 MUONRecoCheck.C:157
 MUONRecoCheck.C:158
 MUONRecoCheck.C:159
 MUONRecoCheck.C:160
 MUONRecoCheck.C:161
 MUONRecoCheck.C:162
 MUONRecoCheck.C:163
 MUONRecoCheck.C:164
 MUONRecoCheck.C:165
 MUONRecoCheck.C:166
 MUONRecoCheck.C:167
 MUONRecoCheck.C:168
 MUONRecoCheck.C:169
 MUONRecoCheck.C:170
 MUONRecoCheck.C:171
 MUONRecoCheck.C:172
 MUONRecoCheck.C:173
 MUONRecoCheck.C:174
 MUONRecoCheck.C:175
 MUONRecoCheck.C:176
 MUONRecoCheck.C:177
 MUONRecoCheck.C:178
 MUONRecoCheck.C:179
 MUONRecoCheck.C:180
 MUONRecoCheck.C:181
 MUONRecoCheck.C:182
 MUONRecoCheck.C:183
 MUONRecoCheck.C:184
 MUONRecoCheck.C:185
 MUONRecoCheck.C:186
 MUONRecoCheck.C:187
 MUONRecoCheck.C:188
 MUONRecoCheck.C:189
 MUONRecoCheck.C:190
 MUONRecoCheck.C:191
 MUONRecoCheck.C:192
 MUONRecoCheck.C:193
 MUONRecoCheck.C:194
 MUONRecoCheck.C:195
 MUONRecoCheck.C:196
 MUONRecoCheck.C:197
 MUONRecoCheck.C:198
 MUONRecoCheck.C:199
 MUONRecoCheck.C:200
 MUONRecoCheck.C:201
 MUONRecoCheck.C:202
 MUONRecoCheck.C:203
 MUONRecoCheck.C:204
 MUONRecoCheck.C:205
 MUONRecoCheck.C:206
 MUONRecoCheck.C:207
 MUONRecoCheck.C:208
 MUONRecoCheck.C:209
 MUONRecoCheck.C:210
 MUONRecoCheck.C:211
 MUONRecoCheck.C:212
 MUONRecoCheck.C:213
 MUONRecoCheck.C:214
 MUONRecoCheck.C:215
 MUONRecoCheck.C:216
 MUONRecoCheck.C:217
 MUONRecoCheck.C:218
 MUONRecoCheck.C:219
 MUONRecoCheck.C:220
 MUONRecoCheck.C:221
 MUONRecoCheck.C:222
 MUONRecoCheck.C:223
 MUONRecoCheck.C:224
 MUONRecoCheck.C:225
 MUONRecoCheck.C:226
 MUONRecoCheck.C:227
 MUONRecoCheck.C:228
 MUONRecoCheck.C:229
 MUONRecoCheck.C:230
 MUONRecoCheck.C:231
 MUONRecoCheck.C:232
 MUONRecoCheck.C:233
 MUONRecoCheck.C:234
 MUONRecoCheck.C:235
 MUONRecoCheck.C:236
 MUONRecoCheck.C:237
 MUONRecoCheck.C:238
 MUONRecoCheck.C:239
 MUONRecoCheck.C:240
 MUONRecoCheck.C:241
 MUONRecoCheck.C:242
 MUONRecoCheck.C:243
 MUONRecoCheck.C:244
 MUONRecoCheck.C:245
 MUONRecoCheck.C:246
 MUONRecoCheck.C:247
 MUONRecoCheck.C:248
 MUONRecoCheck.C:249
 MUONRecoCheck.C:250
 MUONRecoCheck.C:251
 MUONRecoCheck.C:252
 MUONRecoCheck.C:253
 MUONRecoCheck.C:254
 MUONRecoCheck.C:255
 MUONRecoCheck.C:256
 MUONRecoCheck.C:257
 MUONRecoCheck.C:258
 MUONRecoCheck.C:259
 MUONRecoCheck.C:260
 MUONRecoCheck.C:261
 MUONRecoCheck.C:262
 MUONRecoCheck.C:263
 MUONRecoCheck.C:264
 MUONRecoCheck.C:265
 MUONRecoCheck.C:266
 MUONRecoCheck.C:267
 MUONRecoCheck.C:268
 MUONRecoCheck.C:269
 MUONRecoCheck.C:270
 MUONRecoCheck.C:271
 MUONRecoCheck.C:272
 MUONRecoCheck.C:273
 MUONRecoCheck.C:274
 MUONRecoCheck.C:275
 MUONRecoCheck.C:276
 MUONRecoCheck.C:277
 MUONRecoCheck.C:278
 MUONRecoCheck.C:279
 MUONRecoCheck.C:280
 MUONRecoCheck.C:281
 MUONRecoCheck.C:282
 MUONRecoCheck.C:283
 MUONRecoCheck.C:284
 MUONRecoCheck.C:285
 MUONRecoCheck.C:286
 MUONRecoCheck.C:287
 MUONRecoCheck.C:288
 MUONRecoCheck.C:289
 MUONRecoCheck.C:290
 MUONRecoCheck.C:291
 MUONRecoCheck.C:292
 MUONRecoCheck.C:293
 MUONRecoCheck.C:294
 MUONRecoCheck.C:295
 MUONRecoCheck.C:296
 MUONRecoCheck.C:297
 MUONRecoCheck.C:298
 MUONRecoCheck.C:299
 MUONRecoCheck.C:300
 MUONRecoCheck.C:301
 MUONRecoCheck.C:302
 MUONRecoCheck.C:303
 MUONRecoCheck.C:304
 MUONRecoCheck.C:305
 MUONRecoCheck.C:306
 MUONRecoCheck.C:307
 MUONRecoCheck.C:308
 MUONRecoCheck.C:309
 MUONRecoCheck.C:310
 MUONRecoCheck.C:311
 MUONRecoCheck.C:312
 MUONRecoCheck.C:313
 MUONRecoCheck.C:314
 MUONRecoCheck.C:315
 MUONRecoCheck.C:316
 MUONRecoCheck.C:317
 MUONRecoCheck.C:318
 MUONRecoCheck.C:319
 MUONRecoCheck.C:320
 MUONRecoCheck.C:321
 MUONRecoCheck.C:322
 MUONRecoCheck.C:323
 MUONRecoCheck.C:324
 MUONRecoCheck.C:325
 MUONRecoCheck.C:326
 MUONRecoCheck.C:327
 MUONRecoCheck.C:328
 MUONRecoCheck.C:329
 MUONRecoCheck.C:330
 MUONRecoCheck.C:331
 MUONRecoCheck.C:332
 MUONRecoCheck.C:333
 MUONRecoCheck.C:334
 MUONRecoCheck.C:335
 MUONRecoCheck.C:336
 MUONRecoCheck.C:337
 MUONRecoCheck.C:338
 MUONRecoCheck.C:339
 MUONRecoCheck.C:340
 MUONRecoCheck.C:341
 MUONRecoCheck.C:342
 MUONRecoCheck.C:343
 MUONRecoCheck.C:344
 MUONRecoCheck.C:345
 MUONRecoCheck.C:346
 MUONRecoCheck.C:347
 MUONRecoCheck.C:348
 MUONRecoCheck.C:349
 MUONRecoCheck.C:350
 MUONRecoCheck.C:351
 MUONRecoCheck.C:352
 MUONRecoCheck.C:353
 MUONRecoCheck.C:354
 MUONRecoCheck.C:355
 MUONRecoCheck.C:356
 MUONRecoCheck.C:357
 MUONRecoCheck.C:358
 MUONRecoCheck.C:359
 MUONRecoCheck.C:360
 MUONRecoCheck.C:361
 MUONRecoCheck.C:362
 MUONRecoCheck.C:363
 MUONRecoCheck.C:364
 MUONRecoCheck.C:365
 MUONRecoCheck.C:366
 MUONRecoCheck.C:367
 MUONRecoCheck.C:368
 MUONRecoCheck.C:369
 MUONRecoCheck.C:370
 MUONRecoCheck.C:371
 MUONRecoCheck.C:372
 MUONRecoCheck.C:373
 MUONRecoCheck.C:374
 MUONRecoCheck.C:375
 MUONRecoCheck.C:376
 MUONRecoCheck.C:377
 MUONRecoCheck.C:378
 MUONRecoCheck.C:379
 MUONRecoCheck.C:380
 MUONRecoCheck.C:381
 MUONRecoCheck.C:382
 MUONRecoCheck.C:383
 MUONRecoCheck.C:384
 MUONRecoCheck.C:385
 MUONRecoCheck.C:386
 MUONRecoCheck.C:387
 MUONRecoCheck.C:388
 MUONRecoCheck.C:389
 MUONRecoCheck.C:390
 MUONRecoCheck.C:391
 MUONRecoCheck.C:392
 MUONRecoCheck.C:393
 MUONRecoCheck.C:394
 MUONRecoCheck.C:395
 MUONRecoCheck.C:396
 MUONRecoCheck.C:397
 MUONRecoCheck.C:398
 MUONRecoCheck.C:399
 MUONRecoCheck.C:400
 MUONRecoCheck.C:401
 MUONRecoCheck.C:402
 MUONRecoCheck.C:403
 MUONRecoCheck.C:404
 MUONRecoCheck.C:405
 MUONRecoCheck.C:406
 MUONRecoCheck.C:407
 MUONRecoCheck.C:408
 MUONRecoCheck.C:409
 MUONRecoCheck.C:410
 MUONRecoCheck.C:411
 MUONRecoCheck.C:412
 MUONRecoCheck.C:413
 MUONRecoCheck.C:414
 MUONRecoCheck.C:415
 MUONRecoCheck.C:416
 MUONRecoCheck.C:417
 MUONRecoCheck.C:418
 MUONRecoCheck.C:419
 MUONRecoCheck.C:420
 MUONRecoCheck.C:421
 MUONRecoCheck.C:422
 MUONRecoCheck.C:423
 MUONRecoCheck.C:424
 MUONRecoCheck.C:425
 MUONRecoCheck.C:426
 MUONRecoCheck.C:427
 MUONRecoCheck.C:428
 MUONRecoCheck.C:429
 MUONRecoCheck.C:430
 MUONRecoCheck.C:431
 MUONRecoCheck.C:432
 MUONRecoCheck.C:433
 MUONRecoCheck.C:434
 MUONRecoCheck.C:435
 MUONRecoCheck.C:436
 MUONRecoCheck.C:437
 MUONRecoCheck.C:438
 MUONRecoCheck.C:439
 MUONRecoCheck.C:440
 MUONRecoCheck.C:441
 MUONRecoCheck.C:442
 MUONRecoCheck.C:443
 MUONRecoCheck.C:444
 MUONRecoCheck.C:445
 MUONRecoCheck.C:446
 MUONRecoCheck.C:447
 MUONRecoCheck.C:448
 MUONRecoCheck.C:449
 MUONRecoCheck.C:450
 MUONRecoCheck.C:451
 MUONRecoCheck.C:452
 MUONRecoCheck.C:453
 MUONRecoCheck.C:454
 MUONRecoCheck.C:455
 MUONRecoCheck.C:456
 MUONRecoCheck.C:457
 MUONRecoCheck.C:458
 MUONRecoCheck.C:459
 MUONRecoCheck.C:460
 MUONRecoCheck.C:461
 MUONRecoCheck.C:462
 MUONRecoCheck.C:463
 MUONRecoCheck.C:464
 MUONRecoCheck.C:465
 MUONRecoCheck.C:466
 MUONRecoCheck.C:467
 MUONRecoCheck.C:468
 MUONRecoCheck.C:469
 MUONRecoCheck.C:470
 MUONRecoCheck.C:471
 MUONRecoCheck.C:472
 MUONRecoCheck.C:473
 MUONRecoCheck.C:474
 MUONRecoCheck.C:475
 MUONRecoCheck.C:476
 MUONRecoCheck.C:477
 MUONRecoCheck.C:478
 MUONRecoCheck.C:479
 MUONRecoCheck.C:480
 MUONRecoCheck.C:481
 MUONRecoCheck.C:482
 MUONRecoCheck.C:483
 MUONRecoCheck.C:484
 MUONRecoCheck.C:485
 MUONRecoCheck.C:486
 MUONRecoCheck.C:487
 MUONRecoCheck.C:488
 MUONRecoCheck.C:489
 MUONRecoCheck.C:490
 MUONRecoCheck.C:491
 MUONRecoCheck.C:492
 MUONRecoCheck.C:493
 MUONRecoCheck.C:494
 MUONRecoCheck.C:495
 MUONRecoCheck.C:496
 MUONRecoCheck.C:497
 MUONRecoCheck.C:498
 MUONRecoCheck.C:499
 MUONRecoCheck.C:500
 MUONRecoCheck.C:501
 MUONRecoCheck.C:502
 MUONRecoCheck.C:503
 MUONRecoCheck.C:504
 MUONRecoCheck.C:505
 MUONRecoCheck.C:506
 MUONRecoCheck.C:507
 MUONRecoCheck.C:508
 MUONRecoCheck.C:509
 MUONRecoCheck.C:510
 MUONRecoCheck.C:511
 MUONRecoCheck.C:512
 MUONRecoCheck.C:513
 MUONRecoCheck.C:514
 MUONRecoCheck.C:515
 MUONRecoCheck.C:516
 MUONRecoCheck.C:517
 MUONRecoCheck.C:518
 MUONRecoCheck.C:519
 MUONRecoCheck.C:520
 MUONRecoCheck.C:521
 MUONRecoCheck.C:522
 MUONRecoCheck.C:523
 MUONRecoCheck.C:524
 MUONRecoCheck.C:525
 MUONRecoCheck.C:526
 MUONRecoCheck.C:527
 MUONRecoCheck.C:528
 MUONRecoCheck.C:529
 MUONRecoCheck.C:530
 MUONRecoCheck.C:531
 MUONRecoCheck.C:532
 MUONRecoCheck.C:533
 MUONRecoCheck.C:534
 MUONRecoCheck.C:535
 MUONRecoCheck.C:536
 MUONRecoCheck.C:537
 MUONRecoCheck.C:538
 MUONRecoCheck.C:539
 MUONRecoCheck.C:540
 MUONRecoCheck.C:541
 MUONRecoCheck.C:542
 MUONRecoCheck.C:543
 MUONRecoCheck.C:544
 MUONRecoCheck.C:545
 MUONRecoCheck.C:546
 MUONRecoCheck.C:547
 MUONRecoCheck.C:548
 MUONRecoCheck.C:549
 MUONRecoCheck.C:550
 MUONRecoCheck.C:551
 MUONRecoCheck.C:552
 MUONRecoCheck.C:553
 MUONRecoCheck.C:554
 MUONRecoCheck.C:555
 MUONRecoCheck.C:556
 MUONRecoCheck.C:557
 MUONRecoCheck.C:558
 MUONRecoCheck.C:559
 MUONRecoCheck.C:560
 MUONRecoCheck.C:561
 MUONRecoCheck.C:562
 MUONRecoCheck.C:563
 MUONRecoCheck.C:564
 MUONRecoCheck.C:565
 MUONRecoCheck.C:566
 MUONRecoCheck.C:567
 MUONRecoCheck.C:568
 MUONRecoCheck.C:569
 MUONRecoCheck.C:570
 MUONRecoCheck.C:571
 MUONRecoCheck.C:572
 MUONRecoCheck.C:573
 MUONRecoCheck.C:574
 MUONRecoCheck.C:575
 MUONRecoCheck.C:576
 MUONRecoCheck.C:577
 MUONRecoCheck.C:578
 MUONRecoCheck.C:579
 MUONRecoCheck.C:580
 MUONRecoCheck.C:581
 MUONRecoCheck.C:582
 MUONRecoCheck.C:583
 MUONRecoCheck.C:584
 MUONRecoCheck.C:585
 MUONRecoCheck.C:586
 MUONRecoCheck.C:587
 MUONRecoCheck.C:588
 MUONRecoCheck.C:589
 MUONRecoCheck.C:590
 MUONRecoCheck.C:591
 MUONRecoCheck.C:592
 MUONRecoCheck.C:593
 MUONRecoCheck.C:594
 MUONRecoCheck.C:595
 MUONRecoCheck.C:596
 MUONRecoCheck.C:597
 MUONRecoCheck.C:598
 MUONRecoCheck.C:599
 MUONRecoCheck.C:600
 MUONRecoCheck.C:601
 MUONRecoCheck.C:602
 MUONRecoCheck.C:603
 MUONRecoCheck.C:604
 MUONRecoCheck.C:605
 MUONRecoCheck.C:606
 MUONRecoCheck.C:607
 MUONRecoCheck.C:608
 MUONRecoCheck.C:609
 MUONRecoCheck.C:610
 MUONRecoCheck.C:611
 MUONRecoCheck.C:612
 MUONRecoCheck.C:613
 MUONRecoCheck.C:614
 MUONRecoCheck.C:615
 MUONRecoCheck.C:616
 MUONRecoCheck.C:617
 MUONRecoCheck.C:618
 MUONRecoCheck.C:619
 MUONRecoCheck.C:620
 MUONRecoCheck.C:621
 MUONRecoCheck.C:622
 MUONRecoCheck.C:623
 MUONRecoCheck.C:624
 MUONRecoCheck.C:625
 MUONRecoCheck.C:626
 MUONRecoCheck.C:627
 MUONRecoCheck.C:628
 MUONRecoCheck.C:629
 MUONRecoCheck.C:630
 MUONRecoCheck.C:631
 MUONRecoCheck.C:632
 MUONRecoCheck.C:633
 MUONRecoCheck.C:634
 MUONRecoCheck.C:635
 MUONRecoCheck.C:636
 MUONRecoCheck.C:637
 MUONRecoCheck.C:638
 MUONRecoCheck.C:639
 MUONRecoCheck.C:640
 MUONRecoCheck.C:641
 MUONRecoCheck.C:642
 MUONRecoCheck.C:643
 MUONRecoCheck.C:644
 MUONRecoCheck.C:645
 MUONRecoCheck.C:646
 MUONRecoCheck.C:647
 MUONRecoCheck.C:648
 MUONRecoCheck.C:649
 MUONRecoCheck.C:650
 MUONRecoCheck.C:651
 MUONRecoCheck.C:652
 MUONRecoCheck.C:653
 MUONRecoCheck.C:654
 MUONRecoCheck.C:655
 MUONRecoCheck.C:656
 MUONRecoCheck.C:657
 MUONRecoCheck.C:658
 MUONRecoCheck.C:659
 MUONRecoCheck.C:660
 MUONRecoCheck.C:661
 MUONRecoCheck.C:662
 MUONRecoCheck.C:663
 MUONRecoCheck.C:664
 MUONRecoCheck.C:665
 MUONRecoCheck.C:666
 MUONRecoCheck.C:667
 MUONRecoCheck.C:668
 MUONRecoCheck.C:669
 MUONRecoCheck.C:670
 MUONRecoCheck.C:671
 MUONRecoCheck.C:672
 MUONRecoCheck.C:673
 MUONRecoCheck.C:674
 MUONRecoCheck.C:675
 MUONRecoCheck.C:676
 MUONRecoCheck.C:677
 MUONRecoCheck.C:678
 MUONRecoCheck.C:679
 MUONRecoCheck.C:680
 MUONRecoCheck.C:681
 MUONRecoCheck.C:682
 MUONRecoCheck.C:683
 MUONRecoCheck.C:684
 MUONRecoCheck.C:685
 MUONRecoCheck.C:686
 MUONRecoCheck.C:687
 MUONRecoCheck.C:688
 MUONRecoCheck.C:689
 MUONRecoCheck.C:690
 MUONRecoCheck.C:691
 MUONRecoCheck.C:692
 MUONRecoCheck.C:693
 MUONRecoCheck.C:694
 MUONRecoCheck.C:695
 MUONRecoCheck.C:696
 MUONRecoCheck.C:697
 MUONRecoCheck.C:698
 MUONRecoCheck.C:699
 MUONRecoCheck.C:700
 MUONRecoCheck.C:701
 MUONRecoCheck.C:702
 MUONRecoCheck.C:703
 MUONRecoCheck.C:704
 MUONRecoCheck.C:705
 MUONRecoCheck.C:706
 MUONRecoCheck.C:707
 MUONRecoCheck.C:708
 MUONRecoCheck.C:709
 MUONRecoCheck.C:710
 MUONRecoCheck.C:711
 MUONRecoCheck.C:712
 MUONRecoCheck.C:713
 MUONRecoCheck.C:714
 MUONRecoCheck.C:715
 MUONRecoCheck.C:716
 MUONRecoCheck.C:717
 MUONRecoCheck.C:718
 MUONRecoCheck.C:719
 MUONRecoCheck.C:720
 MUONRecoCheck.C:721
 MUONRecoCheck.C:722
 MUONRecoCheck.C:723
 MUONRecoCheck.C:724
 MUONRecoCheck.C:725
 MUONRecoCheck.C:726
 MUONRecoCheck.C:727
 MUONRecoCheck.C:728
 MUONRecoCheck.C:729
 MUONRecoCheck.C:730
 MUONRecoCheck.C:731
 MUONRecoCheck.C:732
 MUONRecoCheck.C:733
 MUONRecoCheck.C:734
 MUONRecoCheck.C:735
 MUONRecoCheck.C:736
 MUONRecoCheck.C:737
 MUONRecoCheck.C:738
 MUONRecoCheck.C:739
 MUONRecoCheck.C:740
 MUONRecoCheck.C:741
 MUONRecoCheck.C:742
 MUONRecoCheck.C:743
 MUONRecoCheck.C:744
 MUONRecoCheck.C:745
 MUONRecoCheck.C:746
 MUONRecoCheck.C:747
 MUONRecoCheck.C:748
 MUONRecoCheck.C:749
 MUONRecoCheck.C:750
 MUONRecoCheck.C:751
 MUONRecoCheck.C:752
 MUONRecoCheck.C:753
 MUONRecoCheck.C:754
 MUONRecoCheck.C:755
 MUONRecoCheck.C:756
 MUONRecoCheck.C:757
 MUONRecoCheck.C:758
 MUONRecoCheck.C:759
 MUONRecoCheck.C:760
 MUONRecoCheck.C:761
 MUONRecoCheck.C:762
 MUONRecoCheck.C:763
 MUONRecoCheck.C:764
 MUONRecoCheck.C:765
 MUONRecoCheck.C:766
 MUONRecoCheck.C:767
 MUONRecoCheck.C:768
 MUONRecoCheck.C:769
 MUONRecoCheck.C:770
 MUONRecoCheck.C:771
 MUONRecoCheck.C:772
 MUONRecoCheck.C:773
 MUONRecoCheck.C:774
 MUONRecoCheck.C:775
 MUONRecoCheck.C:776
 MUONRecoCheck.C:777
 MUONRecoCheck.C:778
 MUONRecoCheck.C:779
 MUONRecoCheck.C:780
 MUONRecoCheck.C:781
 MUONRecoCheck.C:782
 MUONRecoCheck.C:783
 MUONRecoCheck.C:784
 MUONRecoCheck.C:785
 MUONRecoCheck.C:786
 MUONRecoCheck.C:787
 MUONRecoCheck.C:788
 MUONRecoCheck.C:789
 MUONRecoCheck.C:790
 MUONRecoCheck.C:791
 MUONRecoCheck.C:792
 MUONRecoCheck.C:793
 MUONRecoCheck.C:794
 MUONRecoCheck.C:795
 MUONRecoCheck.C:796
 MUONRecoCheck.C:797
 MUONRecoCheck.C:798
 MUONRecoCheck.C:799
 MUONRecoCheck.C:800
 MUONRecoCheck.C:801
 MUONRecoCheck.C:802
 MUONRecoCheck.C:803
 MUONRecoCheck.C:804
 MUONRecoCheck.C:805
 MUONRecoCheck.C:806
 MUONRecoCheck.C:807
 MUONRecoCheck.C:808
 MUONRecoCheck.C:809
 MUONRecoCheck.C:810
 MUONRecoCheck.C:811
 MUONRecoCheck.C:812
 MUONRecoCheck.C:813
 MUONRecoCheck.C:814
 MUONRecoCheck.C:815
 MUONRecoCheck.C:816
 MUONRecoCheck.C:817
 MUONRecoCheck.C:818
 MUONRecoCheck.C:819
 MUONRecoCheck.C:820
 MUONRecoCheck.C:821
 MUONRecoCheck.C:822
 MUONRecoCheck.C:823
 MUONRecoCheck.C:824
 MUONRecoCheck.C:825
 MUONRecoCheck.C:826
 MUONRecoCheck.C:827
 MUONRecoCheck.C:828
 MUONRecoCheck.C:829
 MUONRecoCheck.C:830
 MUONRecoCheck.C:831
 MUONRecoCheck.C:832
 MUONRecoCheck.C:833
 MUONRecoCheck.C:834
 MUONRecoCheck.C:835
 MUONRecoCheck.C:836
 MUONRecoCheck.C:837
 MUONRecoCheck.C:838
 MUONRecoCheck.C:839
 MUONRecoCheck.C:840
 MUONRecoCheck.C:841
 MUONRecoCheck.C:842
 MUONRecoCheck.C:843
 MUONRecoCheck.C:844
 MUONRecoCheck.C:845
 MUONRecoCheck.C:846
 MUONRecoCheck.C:847
 MUONRecoCheck.C:848
 MUONRecoCheck.C:849
 MUONRecoCheck.C:850
 MUONRecoCheck.C:851
 MUONRecoCheck.C:852
 MUONRecoCheck.C:853
 MUONRecoCheck.C:854
 MUONRecoCheck.C:855
 MUONRecoCheck.C:856
 MUONRecoCheck.C:857
 MUONRecoCheck.C:858
 MUONRecoCheck.C:859
 MUONRecoCheck.C:860
 MUONRecoCheck.C:861
 MUONRecoCheck.C:862
 MUONRecoCheck.C:863
 MUONRecoCheck.C:864
 MUONRecoCheck.C:865
 MUONRecoCheck.C:866
 MUONRecoCheck.C:867
 MUONRecoCheck.C:868
 MUONRecoCheck.C:869
 MUONRecoCheck.C:870
 MUONRecoCheck.C:871
 MUONRecoCheck.C:872
 MUONRecoCheck.C:873
 MUONRecoCheck.C:874
 MUONRecoCheck.C:875
 MUONRecoCheck.C:876
 MUONRecoCheck.C:877
 MUONRecoCheck.C:878
 MUONRecoCheck.C:879
 MUONRecoCheck.C:880
 MUONRecoCheck.C:881
 MUONRecoCheck.C:882
 MUONRecoCheck.C:883
 MUONRecoCheck.C:884
 MUONRecoCheck.C:885
 MUONRecoCheck.C:886
 MUONRecoCheck.C:887
 MUONRecoCheck.C:888
 MUONRecoCheck.C:889
 MUONRecoCheck.C:890
 MUONRecoCheck.C:891
 MUONRecoCheck.C:892
 MUONRecoCheck.C:893
 MUONRecoCheck.C:894
 MUONRecoCheck.C:895
 MUONRecoCheck.C:896
 MUONRecoCheck.C:897
 MUONRecoCheck.C:898
 MUONRecoCheck.C:899
 MUONRecoCheck.C:900
 MUONRecoCheck.C:901
 MUONRecoCheck.C:902
 MUONRecoCheck.C:903
 MUONRecoCheck.C:904
 MUONRecoCheck.C:905
 MUONRecoCheck.C:906
 MUONRecoCheck.C:907
 MUONRecoCheck.C:908
 MUONRecoCheck.C:909
 MUONRecoCheck.C:910
 MUONRecoCheck.C:911
 MUONRecoCheck.C:912
 MUONRecoCheck.C:913
 MUONRecoCheck.C:914
 MUONRecoCheck.C:915
 MUONRecoCheck.C:916
 MUONRecoCheck.C:917
 MUONRecoCheck.C:918
 MUONRecoCheck.C:919
 MUONRecoCheck.C:920
 MUONRecoCheck.C:921
 MUONRecoCheck.C:922
 MUONRecoCheck.C:923
 MUONRecoCheck.C:924
 MUONRecoCheck.C:925
 MUONRecoCheck.C:926
 MUONRecoCheck.C:927
 MUONRecoCheck.C:928
 MUONRecoCheck.C:929
 MUONRecoCheck.C:930
 MUONRecoCheck.C:931
 MUONRecoCheck.C:932
 MUONRecoCheck.C:933
 MUONRecoCheck.C:934
 MUONRecoCheck.C:935
 MUONRecoCheck.C:936
 MUONRecoCheck.C:937
 MUONRecoCheck.C:938
 MUONRecoCheck.C:939
 MUONRecoCheck.C:940
 MUONRecoCheck.C:941
 MUONRecoCheck.C:942
 MUONRecoCheck.C:943
 MUONRecoCheck.C:944
 MUONRecoCheck.C:945
 MUONRecoCheck.C:946
 MUONRecoCheck.C:947
 MUONRecoCheck.C:948
 MUONRecoCheck.C:949
 MUONRecoCheck.C:950
 MUONRecoCheck.C:951
 MUONRecoCheck.C:952
 MUONRecoCheck.C:953
 MUONRecoCheck.C:954
 MUONRecoCheck.C:955
 MUONRecoCheck.C:956
 MUONRecoCheck.C:957
 MUONRecoCheck.C:958
 MUONRecoCheck.C:959
 MUONRecoCheck.C:960
 MUONRecoCheck.C:961
 MUONRecoCheck.C:962
 MUONRecoCheck.C:963
 MUONRecoCheck.C:964
 MUONRecoCheck.C:965
 MUONRecoCheck.C:966
 MUONRecoCheck.C:967
 MUONRecoCheck.C:968
 MUONRecoCheck.C:969
 MUONRecoCheck.C:970
 MUONRecoCheck.C:971
 MUONRecoCheck.C:972
 MUONRecoCheck.C:973
 MUONRecoCheck.C:974
 MUONRecoCheck.C:975
 MUONRecoCheck.C:976
 MUONRecoCheck.C:977
 MUONRecoCheck.C:978
 MUONRecoCheck.C:979
 MUONRecoCheck.C:980
 MUONRecoCheck.C:981
 MUONRecoCheck.C:982
 MUONRecoCheck.C:983
 MUONRecoCheck.C:984
 MUONRecoCheck.C:985
 MUONRecoCheck.C:986
 MUONRecoCheck.C:987
 MUONRecoCheck.C:988
 MUONRecoCheck.C:989
 MUONRecoCheck.C:990
 MUONRecoCheck.C:991
 MUONRecoCheck.C:992
 MUONRecoCheck.C:993
 MUONRecoCheck.C:994
 MUONRecoCheck.C:995
 MUONRecoCheck.C:996
 MUONRecoCheck.C:997
 MUONRecoCheck.C:998
 MUONRecoCheck.C:999
 MUONRecoCheck.C:1000
 MUONRecoCheck.C:1001
 MUONRecoCheck.C:1002
 MUONRecoCheck.C:1003
 MUONRecoCheck.C:1004
 MUONRecoCheck.C:1005
 MUONRecoCheck.C:1006
 MUONRecoCheck.C:1007
 MUONRecoCheck.C:1008
 MUONRecoCheck.C:1009
 MUONRecoCheck.C:1010
 MUONRecoCheck.C:1011
 MUONRecoCheck.C:1012
 MUONRecoCheck.C:1013
 MUONRecoCheck.C:1014
 MUONRecoCheck.C:1015
 MUONRecoCheck.C:1016
 MUONRecoCheck.C:1017
 MUONRecoCheck.C:1018
 MUONRecoCheck.C:1019
 MUONRecoCheck.C:1020
 MUONRecoCheck.C:1021
 MUONRecoCheck.C:1022
 MUONRecoCheck.C:1023
 MUONRecoCheck.C:1024
 MUONRecoCheck.C:1025
 MUONRecoCheck.C:1026
 MUONRecoCheck.C:1027
 MUONRecoCheck.C:1028
 MUONRecoCheck.C:1029
 MUONRecoCheck.C:1030
 MUONRecoCheck.C:1031
 MUONRecoCheck.C:1032
 MUONRecoCheck.C:1033
 MUONRecoCheck.C:1034
 MUONRecoCheck.C:1035
 MUONRecoCheck.C:1036
 MUONRecoCheck.C:1037
 MUONRecoCheck.C:1038
 MUONRecoCheck.C:1039
 MUONRecoCheck.C:1040
 MUONRecoCheck.C:1041
 MUONRecoCheck.C:1042
 MUONRecoCheck.C:1043
 MUONRecoCheck.C:1044
 MUONRecoCheck.C:1045
 MUONRecoCheck.C:1046
 MUONRecoCheck.C:1047
 MUONRecoCheck.C:1048
 MUONRecoCheck.C:1049
 MUONRecoCheck.C:1050
 MUONRecoCheck.C:1051
 MUONRecoCheck.C:1052
 MUONRecoCheck.C:1053
 MUONRecoCheck.C:1054
 MUONRecoCheck.C:1055
 MUONRecoCheck.C:1056
 MUONRecoCheck.C:1057
 MUONRecoCheck.C:1058
 MUONRecoCheck.C:1059
 MUONRecoCheck.C:1060
 MUONRecoCheck.C:1061
 MUONRecoCheck.C:1062
 MUONRecoCheck.C:1063
 MUONRecoCheck.C:1064
 MUONRecoCheck.C:1065
 MUONRecoCheck.C:1066
 MUONRecoCheck.C:1067
 MUONRecoCheck.C:1068
 MUONRecoCheck.C:1069
 MUONRecoCheck.C:1070
 MUONRecoCheck.C:1071
 MUONRecoCheck.C:1072
 MUONRecoCheck.C:1073
 MUONRecoCheck.C:1074
 MUONRecoCheck.C:1075
 MUONRecoCheck.C:1076
 MUONRecoCheck.C:1077
 MUONRecoCheck.C:1078
 MUONRecoCheck.C:1079
 MUONRecoCheck.C:1080
 MUONRecoCheck.C:1081
 MUONRecoCheck.C:1082
 MUONRecoCheck.C:1083
 MUONRecoCheck.C:1084
 MUONRecoCheck.C:1085
 MUONRecoCheck.C:1086
 MUONRecoCheck.C:1087
 MUONRecoCheck.C:1088
 MUONRecoCheck.C:1089
 MUONRecoCheck.C:1090
 MUONRecoCheck.C:1091
 MUONRecoCheck.C:1092
 MUONRecoCheck.C:1093
 MUONRecoCheck.C:1094
 MUONRecoCheck.C:1095
 MUONRecoCheck.C:1096
 MUONRecoCheck.C:1097
 MUONRecoCheck.C:1098
 MUONRecoCheck.C:1099
 MUONRecoCheck.C:1100
 MUONRecoCheck.C:1101
 MUONRecoCheck.C:1102
 MUONRecoCheck.C:1103
 MUONRecoCheck.C:1104
 MUONRecoCheck.C:1105
 MUONRecoCheck.C:1106
 MUONRecoCheck.C:1107
 MUONRecoCheck.C:1108
 MUONRecoCheck.C:1109
 MUONRecoCheck.C:1110
 MUONRecoCheck.C:1111
 MUONRecoCheck.C:1112
 MUONRecoCheck.C:1113
 MUONRecoCheck.C:1114
 MUONRecoCheck.C:1115
 MUONRecoCheck.C:1116
 MUONRecoCheck.C:1117
 MUONRecoCheck.C:1118
 MUONRecoCheck.C:1119
 MUONRecoCheck.C:1120
 MUONRecoCheck.C:1121
 MUONRecoCheck.C:1122
 MUONRecoCheck.C:1123
 MUONRecoCheck.C:1124
 MUONRecoCheck.C:1125
 MUONRecoCheck.C:1126
 MUONRecoCheck.C:1127
 MUONRecoCheck.C:1128
 MUONRecoCheck.C:1129
 MUONRecoCheck.C:1130
 MUONRecoCheck.C:1131
 MUONRecoCheck.C:1132
 MUONRecoCheck.C:1133
 MUONRecoCheck.C:1134
 MUONRecoCheck.C:1135
 MUONRecoCheck.C:1136
 MUONRecoCheck.C:1137
 MUONRecoCheck.C:1138
 MUONRecoCheck.C:1139
 MUONRecoCheck.C:1140
 MUONRecoCheck.C:1141
 MUONRecoCheck.C:1142
 MUONRecoCheck.C:1143
 MUONRecoCheck.C:1144
 MUONRecoCheck.C:1145
 MUONRecoCheck.C:1146
 MUONRecoCheck.C:1147
 MUONRecoCheck.C:1148
 MUONRecoCheck.C:1149
 MUONRecoCheck.C:1150
 MUONRecoCheck.C:1151
 MUONRecoCheck.C:1152
 MUONRecoCheck.C:1153
 MUONRecoCheck.C:1154
 MUONRecoCheck.C:1155
 MUONRecoCheck.C:1156
 MUONRecoCheck.C:1157
 MUONRecoCheck.C:1158
 MUONRecoCheck.C:1159
 MUONRecoCheck.C:1160
 MUONRecoCheck.C:1161
 MUONRecoCheck.C:1162
 MUONRecoCheck.C:1163
 MUONRecoCheck.C:1164
 MUONRecoCheck.C:1165
 MUONRecoCheck.C:1166
 MUONRecoCheck.C:1167
 MUONRecoCheck.C:1168
 MUONRecoCheck.C:1169
 MUONRecoCheck.C:1170
 MUONRecoCheck.C:1171
 MUONRecoCheck.C:1172
 MUONRecoCheck.C:1173
 MUONRecoCheck.C:1174
 MUONRecoCheck.C:1175
 MUONRecoCheck.C:1176
 MUONRecoCheck.C:1177
 MUONRecoCheck.C:1178
 MUONRecoCheck.C:1179
 MUONRecoCheck.C:1180
 MUONRecoCheck.C:1181
 MUONRecoCheck.C:1182
 MUONRecoCheck.C:1183
 MUONRecoCheck.C:1184
 MUONRecoCheck.C:1185
 MUONRecoCheck.C:1186
 MUONRecoCheck.C:1187
 MUONRecoCheck.C:1188
 MUONRecoCheck.C:1189
 MUONRecoCheck.C:1190
 MUONRecoCheck.C:1191
 MUONRecoCheck.C:1192
 MUONRecoCheck.C:1193
 MUONRecoCheck.C:1194
 MUONRecoCheck.C:1195
 MUONRecoCheck.C:1196
 MUONRecoCheck.C:1197
 MUONRecoCheck.C:1198
 MUONRecoCheck.C:1199
 MUONRecoCheck.C:1200
 MUONRecoCheck.C:1201
 MUONRecoCheck.C:1202
 MUONRecoCheck.C:1203
 MUONRecoCheck.C:1204
 MUONRecoCheck.C:1205
 MUONRecoCheck.C:1206
 MUONRecoCheck.C:1207
 MUONRecoCheck.C:1208
 MUONRecoCheck.C:1209
 MUONRecoCheck.C:1210
 MUONRecoCheck.C:1211
 MUONRecoCheck.C:1212
 MUONRecoCheck.C:1213
 MUONRecoCheck.C:1214
 MUONRecoCheck.C:1215
 MUONRecoCheck.C:1216
 MUONRecoCheck.C:1217
 MUONRecoCheck.C:1218
 MUONRecoCheck.C:1219
 MUONRecoCheck.C:1220
 MUONRecoCheck.C:1221
 MUONRecoCheck.C:1222
 MUONRecoCheck.C:1223
 MUONRecoCheck.C:1224
 MUONRecoCheck.C:1225
 MUONRecoCheck.C:1226
 MUONRecoCheck.C:1227
 MUONRecoCheck.C:1228
 MUONRecoCheck.C:1229
 MUONRecoCheck.C:1230
 MUONRecoCheck.C:1231
 MUONRecoCheck.C:1232
 MUONRecoCheck.C:1233
 MUONRecoCheck.C:1234
 MUONRecoCheck.C:1235
 MUONRecoCheck.C:1236
 MUONRecoCheck.C:1237
 MUONRecoCheck.C:1238
 MUONRecoCheck.C:1239
 MUONRecoCheck.C:1240
 MUONRecoCheck.C:1241
 MUONRecoCheck.C:1242
 MUONRecoCheck.C:1243
 MUONRecoCheck.C:1244
 MUONRecoCheck.C:1245
 MUONRecoCheck.C:1246
 MUONRecoCheck.C:1247
 MUONRecoCheck.C:1248
 MUONRecoCheck.C:1249
 MUONRecoCheck.C:1250
 MUONRecoCheck.C:1251
 MUONRecoCheck.C:1252