ROOT logo
#ifndef ALIMILLEPEDE2_H
#define ALIMILLEPEDE2_H

/**********************************************************************************************/
/* General class for alignment with large number of degrees of freedom                        */
/* Based on the original milliped2 by Volker Blobel                                           */
/* http://www.desy.de/~blobel/mptalks.html                                                    */
/*                                                                                            */ 
/* Author: ruben.shahoyan@cern.ch                                                             */
/*                                                                                            */ 
/**********************************************************************************************/

#include <TObject.h>
#include <TString.h>
#include <TTree.h>
#include "AliMinResSolve.h"
#include "AliMillePedeRecord.h"
class TFile;
class AliMatrixSq;
class AliSymMatrix;
class AliRectMatrix;
class AliMatrixSparse;
class AliLog;
class TStopwatch;
class TArrayL;
class TArrayF;


class AliMillePede2: public TObject
{
 public:
  //
  enum {kFailed,kInvert,kNoInversion};    // used global matrix solution methods
  enum {kFixParID=-1};                    // dummy id for fixed param
  //
  AliMillePede2();
  AliMillePede2(const AliMillePede2& src);
  virtual ~AliMillePede2();
  AliMillePede2& operator=(const AliMillePede2& )             {printf("Dummy\n"); return *this;}
  //
  Int_t                InitMille(int nGlo, int nLoc, Int_t lNStdDev=-1,double lResCut=-1., double lResCutInit=-1., const Int_t* regroup=0);
  //
  Int_t                GetNGloPar()                     const {return fNGloPar;}
  Int_t                GetNGloParIni()                  const {return fNGloParIni;}
  const Int_t*         GetRegrouping()                  const {return fkReGroup;}
  Int_t                GetNLocPar()                     const {return fNLocPar;}
  Long_t               GetNLocalEquations()             const {return fNLocEquations;}
  Int_t                GetCurrentIteration()            const {return fIter;}
  Int_t                GetNMaxIterations()              const {return fMaxIter;}
  Int_t                GetNStdDev()                     const {return fNStdDev;} 
  Int_t                GetNGlobalConstraints()          const {return fNGloConstraints;}
  Int_t                GetNLagrangeConstraints()        const {return fNLagrangeConstraints;}
  Long_t               GetNLocalFits()                  const {return fNLocFits;}
  Long_t               GetNLocalFitsRejected()          const {return fNLocFitsRejected;}
  Int_t                GetNGlobalsFixed()               const {return fNGloFix;}
  Int_t                GetGlobalSolveStatus()           const {return fGloSolveStatus;}
  Float_t              GetChi2CutFactor()               const {return fChi2CutFactor;}
  Float_t              GetChi2CutRef()                  const {return fChi2CutRef;}
  Float_t              GetResCurInit()                  const {return fResCutInit;}
  Float_t              GetResCut()                      const {return fResCut;}
  Int_t                GetMinPntValid()                 const {return fMinPntValid;}
  Int_t                GetRGId(Int_t i)                 const {return fkReGroup ? (fkReGroup[i]<0 ? -1:fkReGroup[i]) : i;}       
  Int_t                GetProcessedPoints(Int_t i)      const {int ir=GetRGId(i); return ir<=0 ? 0:fProcPnt[ir];}
  Int_t*               GetProcessedPoints()             const {return fProcPnt;}
  Int_t                GetParamGrID(Int_t i)            const {int ir=GetRGId(i); return ir<=0 ? 0:fParamGrID[ir];}
  //
  AliMatrixSq*         GetGlobalMatrix()                const {return fMatCGlo;}
  AliSymMatrix*        GetLocalMatrix()                 const {return fMatCLoc;}
  Double_t*            GetGlobals()                     const {return fVecBGlo;}
  Double_t*            GetDeltaPars()                   const {return fDeltaPar;}
  Double_t*            GetInitPars()                    const {return fInitPar;}
  Double_t*            GetSigmaPars()                   const {return fSigmaPar;}
  Bool_t*              GetIsLinear()                    const {return fIsLinear;}
  Double_t             GetFinalParam(int i)             const {int ir=GetRGId(i); return ir<0 ? 0:fDeltaPar[ir]+fInitPar[ir];}
  Double_t             GetFinalError(int i)             const {return GetParError(i);}
  Double_t             GetPull(int i)                   const;
  //
  Double_t             GetGlobal(Int_t i)               const {int ir=GetRGId(i); return ir<0 ? 0:fVecBGlo[ir];}
  Double_t             GetInitPar(Int_t i)              const {int ir=GetRGId(i); return ir<0 ? 0:fInitPar[ir];}
  Double_t             GetSigmaPar(Int_t i)             const {int ir=GetRGId(i); return ir<0 ? 0:fSigmaPar[ir];}
  Bool_t               GetIsLinear(Int_t i)             const {int ir=GetRGId(i); return ir<0 ? 0:fIsLinear[ir];}
  static Bool_t        IsGlobalMatSparse()                    {return fgIsMatGloSparse;}
  static Bool_t        IsWeightSigma()                        {return fgWeightSigma;}
  void                 SetWghScale(Double_t wOdd=1,Double_t wEven=1)    {fWghScl[0] = wOdd; fWghScl[1] = wEven;}
  void                 SetUseRecordWeight(Bool_t v=kTRUE)     {fUseRecordWeight=v;}
  Bool_t               GetUseRecordWeight()             const {return fUseRecordWeight;}
  void                 SetMinRecordLength(Int_t v=1)          {fMinRecordLength = v;}
  Int_t                GetMinRecordLength()             const {return fMinRecordLength;}
  //
  void                 SetParamGrID(Int_t grID,Int_t i)       {int ir=GetRGId(i); if(ir<0) return; fParamGrID[ir] = grID; if(fNGroupsSet<grID)fNGroupsSet=grID;}
  void                 SetNGloPar(Int_t n)                    {fNGloPar = n;}
  void                 SetNLocPar(Int_t n)                    {fNLocPar = n;}
  void                 SetNMaxIterations(Int_t n=10)          {fMaxIter = n;}
  void                 SetNStdDev(Int_t n)                    {fNStdDev = n;}
  void                 SetChi2CutFactor(Float_t v)            {fChi2CutFactor = v;}
  void                 SetChi2CutRef(Float_t v)               {fChi2CutRef = v;}
  void                 SetResCurInit(Float_t v)               {fResCutInit = v;}
  void                 SetResCut(Float_t v)                   {fResCut = v;}
  void                 SetMinPntValid(Int_t n)                {fMinPntValid = n>0 ? n:1;}
  static void          SetGlobalMatSparse(Bool_t v=kTRUE)     {fgIsMatGloSparse = v;}
  static void          SetWeightSigma(Bool_t v=kTRUE)         {fgWeightSigma = v;}
  //
  void                 SetInitPars(const Double_t* par);
  void                 SetSigmaPars(const Double_t* par);
  void                 SetInitPar(Int_t i,Double_t par);
  void                 SetSigmaPar(Int_t i,Double_t par);
  //
  Int_t                GlobalFit(Double_t *par=0, Double_t *error=0, Double_t *pull=0);
  Int_t                GlobalFitIteration();
  Int_t                SolveGlobalMatEq();
  static void          SetInvChol(Bool_t v=kTRUE)             {fgInvChol = v;}
  static void          SetMinResPrecondType(Int_t tp=0)       {fgMinResCondType = tp;}
  static void          SetMinResTol(Double_t val=1e-12)       {fgMinResTol = val;}
  static void          SetMinResMaxIter(Int_t val=2000)       {fgMinResMaxIter = val;}
  static void          SetIterSolverType(Int_t val=AliMinResSolve::kSolMinRes) {fgIterSol = val;}
  static void          SetNKrylovV(Int_t val=60)              {fgNKrylovV = val;}
  //
  static Bool_t        GetInvChol()                           {return fgInvChol;}
  static Int_t         GetMinResPrecondType()                 {return fgMinResCondType;}
  static Double_t      GetMinResTol()                         {return fgMinResTol;}
  static Int_t         GetMinResMaxIter()                     {return fgMinResMaxIter;}
  static Int_t         GetIterSolverType()                    {return fgIterSol;}
  static Int_t         GetNKrylovV()                          {return fgNKrylovV;}
  //
  Double_t             GetParError(int iPar)           const;
  Int_t                PrintGlobalParameters()         const;
  void                 SetRejRunList(const UInt_t *runs, Int_t nruns);
  void                 SetAccRunList(const UInt_t *runs, Int_t nruns, const Float_t* wghList=0);
  Bool_t               IsRecordAcceptable();
  //
  //
  Int_t                SetIterations(double lChi2CutFac);

  //
  // constraints
  void                 SetGlobalConstraint(const double *dergb, double val, double sigma=0);
  void                 SetGlobalConstraint(const int *indgb, const double *dergb, int ngb, double val, double sigma=0);
  //
  // processing of the local measurement
  void                 SetRecordRun(Int_t run);
  void                 SetRecordWeight(double wgh);
  void                 SetLocalEquation(double *dergb, double *derlc, double lMeas, double lSigma);
  void                 SetLocalEquation(int *indgb, double *dergb, int ngb, int *indlc, 
					double *derlc,int nlc,double lMeas,double lSigma);
  //
  // manipilation with processed data and costraints records and its buffer
  void                 SetDataRecFName(const char* flname)   {fDataRecFName = flname;}
  const Char_t*        GetDataRecFName()               const {return fDataRecFName.Data();}
  void                 SetConsRecFName(const char* flname)   {fConstrRecFName = flname;}
  const Char_t*        GetConsRecFName()               const {return fConstrRecFName.Data();}
  //
  void   SetRecDataTreeName(const char* name=0)     {fRecDataTreeName = name;   if (fRecDataTreeName.IsNull()) fRecDataTreeName = "AliMillePedeRecords_Data";}
  void   SetRecConsTreeName(const char* name=0)     {fRecConsTreeName = name;   if (fRecConsTreeName.IsNull()) fRecConsTreeName = "AliMillePedeRecords_Consaints";}
  void   SetRecDataBranchName(const char* name=0)   {fRecDataBranchName = name; if (fRecDataBranchName.IsNull()) fRecDataBranchName = "Record_Data";}
  void   SetRecConsBranchName(const char* name=0)   {fRecConsBranchName = name; if (fRecConsBranchName.IsNull()) fRecConsBranchName = "Record_Consaints";}
  const char* GetRecDataTreeName()     const {return fRecDataTreeName.Data();}
  const char* GetRecConsTreeName()     const {return fRecConsTreeName.Data();}
  const char* GetRecDataBranchName()   const {return fRecDataBranchName.Data();}
  const char* GetRecConsBranchName()   const {return fRecConsBranchName.Data();}
  //
  Bool_t               InitDataRecStorage(Bool_t read=kFALSE);
  Bool_t               InitConsRecStorage(Bool_t read=kFALSE);
  Bool_t               ImposeDataRecFile(const char* fname);
  Bool_t               ImposeConsRecFile(const char* fname);
  void                 CloseDataRecStorage();
  void                 CloseConsRecStorage();
  void                 ReadRecordData(Long_t recID);
  void                 ReadRecordConstraint(Long_t recID);
  Bool_t               ReadNextRecordData();
  Bool_t               ReadNextRecordConstraint();
  void                 SaveRecordData();
  void                 SaveRecordConstraint();
  AliMillePedeRecord*  GetRecord()                      const {return fRecord;}
  Long_t               GetSelFirst()                    const {return fSelFirst;}
  Long_t               GetSelLast()                     const {return fSelLast;}
  void                 SetSelFirst(Long_t v)                  {fSelFirst = v;}
  void                 SetSelLast(Long_t v)                   {fSelLast = v;}
  //
  Float_t              Chi2DoFLim(int nSig, int nDoF)   const;
  //
  // aliases for compatibility with millipede1
  void                 SetParSigma(Int_t i,Double_t par)      {SetSigmaPar(i,par);}
  void                 SetGlobalParameters(Double_t *par)     {SetInitPars(par);}
  void                 SetNonLinear(int index, Bool_t v=kTRUE) {int id = GetRGId(index); if (id<0) return; fIsLinear[id] = !v;}
  //
 protected:
  //
  Int_t                LocalFit(double *localParams=0);
  Bool_t               IsZero(Double_t v,Double_t eps=1e-16)   const {return TMath::Abs(v)<eps;}                  
  //
 protected:
  //
  Int_t                 fNLocPar;                        // number of local parameters
  Int_t                 fNGloPar;                        // number of global parameters
  Int_t                 fNGloParIni;                     // number of global parameters before grouping
  Int_t                 fNGloSize;                       // final size of the global matrix (NGloPar+NConstraints)
  //
  Long_t                fNLocEquations;                  // Number of local equations 
  Int_t                 fIter;                           // Current iteration
  Int_t                 fMaxIter;                        // Maximum number of iterations
  Int_t                 fNStdDev;                        // Number of standard deviations for chi2 cut 
  Int_t                 fNGloConstraints;                // Number of constraint equations
  Int_t                 fNLagrangeConstraints;           // Number of constraint equations requiring Lagrange multiplier
  Long_t                fNLocFits;                       // Number of local fits
  Long_t                fNLocFitsRejected;               // Number of local fits rejected
  Int_t                 fNGloFix;                        // Number of globals fixed by user
  Int_t                 fGloSolveStatus;                 // Status of global solver at current step
  //
  Float_t               fChi2CutFactor;                  // Cut factor for chi2 cut to accept local fit 
  Float_t               fChi2CutRef;                     // Reference cut for chi2 cut to accept local fit 
  Float_t               fResCutInit;                     // Cut in residual for first iterartion
  Float_t               fResCut;                         // Cut in residual for other iterartiona
  Int_t                 fMinPntValid;                    // min number of points for global to vary
  //
  Int_t                 fNGroupsSet;                     // number of groups set
  Int_t                *fParamGrID;                      //[fNGloPar] group id for the every parameter
  Int_t                *fProcPnt;                        //[fNGloPar] N of processed points per global variable
  Double_t             *fVecBLoc;                        //[fNLocPar] Vector B local (parameters) 
  Double_t             *fDiagCGlo;                       //[fNGloPar] Initial diagonal elements of C global matrix
  Double_t             *fVecBGlo;                        //! Vector B global (parameters)
  //
  Double_t             *fInitPar;                        //[fNGloPar] Initial global parameters
  Double_t             *fDeltaPar;                       //[fNGloPar] Variation of global parameters
  Double_t             *fSigmaPar;                       //[fNGloPar] Sigma of allowed variation of global parameter
  //
  Bool_t               *fIsLinear;                       //[fNGloPar] Flag for linear parameters
  Bool_t               *fConstrUsed;                     //! Flag for used constraints
  //
  Int_t                *fGlo2CGlo;                       //[fNGloPar] global ID to compressed ID buffer
  Int_t                *fCGlo2Glo;                       //[fNGloPar] compressed ID to global ID buffer
  //
  // Matrices
  AliSymMatrix         *fMatCLoc;                        // Matrix C local
  AliMatrixSq          *fMatCGlo;                        // Matrix C global
  AliRectMatrix        *fMatCGloLoc;                     // Rectangular matrix C g*l 
  Int_t                *fFillIndex;                      //[fNGloPar] auxilary index array for fast matrix fill
  Double_t             *fFillValue;                      //[fNGloPar] auxilary value array for fast matrix fill
  //
  // processed data record bufferization   
  TString               fRecDataTreeName;                // Name of data records tree
  TString               fRecConsTreeName;                // Name of constraints records tree
  TString               fRecDataBranchName;              // Name of data records branch name
  TString               fRecConsBranchName;              // Name of constraints records branch name
  
  TString               fDataRecFName;                   // Name of File for data records               
  AliMillePedeRecord   *fRecord;                         // Buffer of measurements records
  TFile                *fDataRecFile;                    // File of processed measurements records
  TTree                *fTreeData;                       // Tree of processed measurements records
  Int_t                 fRecFileStatus;                  // state of the record file (0-no, 1-read, 2-rw)
  //
  TString               fConstrRecFName;                 // Name of File for constraints records               
  TTree                *fTreeConstr;                     //! Tree of constraint records
  TFile                *fConsRecFile;                    //! File of processed constraints records
  Long_t                fCurrRecDataID;                  // ID of the current data record
  Long_t                fCurrRecConstrID;                // ID of the current constraint record
  Bool_t                fLocFitAdd;                      // Add contribution of carrent track (and not eliminate it)
  Bool_t                fUseRecordWeight;                // force or ignore the record weight
  Int_t                 fMinRecordLength;                // ignore shorter records
  Int_t                 fSelFirst;                       // event selection start
  Int_t                 fSelLast;                        // event selection end
  TArrayL*              fRejRunList;                     // list of runs to reject (if any)
  TArrayL*              fAccRunList;                     // list of runs to select (if any)
  TArrayF*              fAccRunListWgh;                  // optional weights for data of accepted runs (if any)
  Double_t              fRunWgh;                         // run weight
  Double_t              fWghScl[2];                      // optional rescaling for odd/even residual weights (see its usage in LocalFit)
  const Int_t*          fkReGroup;                       // optional regrouping of parameters wrt ID's from the records
  //
  static Bool_t         fgInvChol;                       // Invert global matrix in Cholesky solver
  static Bool_t         fgWeightSigma;                   // weight parameter constraint by statistics
  static Bool_t         fgIsMatGloSparse;                // Type of the global matrix (sparse ...)
  static Int_t          fgMinResCondType;                // Type of the preconditioner for MinRes method 
  static Double_t       fgMinResTol;                     // Tolerance for MinRes solution
  static Int_t          fgMinResMaxIter;                 // Max number of iterations for the MinRes method
  static Int_t          fgIterSol;                       // type of iterative solution: MinRes or FGMRES
  static Int_t          fgNKrylovV;                      // size of Krylov vectors buffer in FGMRES
  //
  ClassDef(AliMillePede2,1)
};

//_____________________________________________________________________________________________
inline void AliMillePede2::ReadRecordData(Long_t recID)       {fTreeData->GetEntry(recID); fCurrRecDataID=recID;}

//_____________________________________________________________________________________________
inline void AliMillePede2::ReadRecordConstraint(Long_t recID) {fTreeConstr->GetEntry(recID); fCurrRecConstrID=recID;}

//_____________________________________________________________________________________________
inline void AliMillePede2::SaveRecordData()                   {fTreeData->Fill(); fRecord->Reset(); fCurrRecDataID++;}

//_____________________________________________________________________________________________
inline void AliMillePede2::SaveRecordConstraint()             {fTreeConstr->Fill(); fRecord->Reset();fCurrRecConstrID++;}

#endif
 AliMillePede2.h:1
 AliMillePede2.h:2
 AliMillePede2.h:3
 AliMillePede2.h:4
 AliMillePede2.h:5
 AliMillePede2.h:6
 AliMillePede2.h:7
 AliMillePede2.h:8
 AliMillePede2.h:9
 AliMillePede2.h:10
 AliMillePede2.h:11
 AliMillePede2.h:12
 AliMillePede2.h:13
 AliMillePede2.h:14
 AliMillePede2.h:15
 AliMillePede2.h:16
 AliMillePede2.h:17
 AliMillePede2.h:18
 AliMillePede2.h:19
 AliMillePede2.h:20
 AliMillePede2.h:21
 AliMillePede2.h:22
 AliMillePede2.h:23
 AliMillePede2.h:24
 AliMillePede2.h:25
 AliMillePede2.h:26
 AliMillePede2.h:27
 AliMillePede2.h:28
 AliMillePede2.h:29
 AliMillePede2.h:30
 AliMillePede2.h:31
 AliMillePede2.h:32
 AliMillePede2.h:33
 AliMillePede2.h:34
 AliMillePede2.h:35
 AliMillePede2.h:36
 AliMillePede2.h:37
 AliMillePede2.h:38
 AliMillePede2.h:39
 AliMillePede2.h:40
 AliMillePede2.h:41
 AliMillePede2.h:42
 AliMillePede2.h:43
 AliMillePede2.h:44
 AliMillePede2.h:45
 AliMillePede2.h:46
 AliMillePede2.h:47
 AliMillePede2.h:48
 AliMillePede2.h:49
 AliMillePede2.h:50
 AliMillePede2.h:51
 AliMillePede2.h:52
 AliMillePede2.h:53
 AliMillePede2.h:54
 AliMillePede2.h:55
 AliMillePede2.h:56
 AliMillePede2.h:57
 AliMillePede2.h:58
 AliMillePede2.h:59
 AliMillePede2.h:60
 AliMillePede2.h:61
 AliMillePede2.h:62
 AliMillePede2.h:63
 AliMillePede2.h:64
 AliMillePede2.h:65
 AliMillePede2.h:66
 AliMillePede2.h:67
 AliMillePede2.h:68
 AliMillePede2.h:69
 AliMillePede2.h:70
 AliMillePede2.h:71
 AliMillePede2.h:72
 AliMillePede2.h:73
 AliMillePede2.h:74
 AliMillePede2.h:75
 AliMillePede2.h:76
 AliMillePede2.h:77
 AliMillePede2.h:78
 AliMillePede2.h:79
 AliMillePede2.h:80
 AliMillePede2.h:81
 AliMillePede2.h:82
 AliMillePede2.h:83
 AliMillePede2.h:84
 AliMillePede2.h:85
 AliMillePede2.h:86
 AliMillePede2.h:87
 AliMillePede2.h:88
 AliMillePede2.h:89
 AliMillePede2.h:90
 AliMillePede2.h:91
 AliMillePede2.h:92
 AliMillePede2.h:93
 AliMillePede2.h:94
 AliMillePede2.h:95
 AliMillePede2.h:96
 AliMillePede2.h:97
 AliMillePede2.h:98
 AliMillePede2.h:99
 AliMillePede2.h:100
 AliMillePede2.h:101
 AliMillePede2.h:102
 AliMillePede2.h:103
 AliMillePede2.h:104
 AliMillePede2.h:105
 AliMillePede2.h:106
 AliMillePede2.h:107
 AliMillePede2.h:108
 AliMillePede2.h:109
 AliMillePede2.h:110
 AliMillePede2.h:111
 AliMillePede2.h:112
 AliMillePede2.h:113
 AliMillePede2.h:114
 AliMillePede2.h:115
 AliMillePede2.h:116
 AliMillePede2.h:117
 AliMillePede2.h:118
 AliMillePede2.h:119
 AliMillePede2.h:120
 AliMillePede2.h:121
 AliMillePede2.h:122
 AliMillePede2.h:123
 AliMillePede2.h:124
 AliMillePede2.h:125
 AliMillePede2.h:126
 AliMillePede2.h:127
 AliMillePede2.h:128
 AliMillePede2.h:129
 AliMillePede2.h:130
 AliMillePede2.h:131
 AliMillePede2.h:132
 AliMillePede2.h:133
 AliMillePede2.h:134
 AliMillePede2.h:135
 AliMillePede2.h:136
 AliMillePede2.h:137
 AliMillePede2.h:138
 AliMillePede2.h:139
 AliMillePede2.h:140
 AliMillePede2.h:141
 AliMillePede2.h:142
 AliMillePede2.h:143
 AliMillePede2.h:144
 AliMillePede2.h:145
 AliMillePede2.h:146
 AliMillePede2.h:147
 AliMillePede2.h:148
 AliMillePede2.h:149
 AliMillePede2.h:150
 AliMillePede2.h:151
 AliMillePede2.h:152
 AliMillePede2.h:153
 AliMillePede2.h:154
 AliMillePede2.h:155
 AliMillePede2.h:156
 AliMillePede2.h:157
 AliMillePede2.h:158
 AliMillePede2.h:159
 AliMillePede2.h:160
 AliMillePede2.h:161
 AliMillePede2.h:162
 AliMillePede2.h:163
 AliMillePede2.h:164
 AliMillePede2.h:165
 AliMillePede2.h:166
 AliMillePede2.h:167
 AliMillePede2.h:168
 AliMillePede2.h:169
 AliMillePede2.h:170
 AliMillePede2.h:171
 AliMillePede2.h:172
 AliMillePede2.h:173
 AliMillePede2.h:174
 AliMillePede2.h:175
 AliMillePede2.h:176
 AliMillePede2.h:177
 AliMillePede2.h:178
 AliMillePede2.h:179
 AliMillePede2.h:180
 AliMillePede2.h:181
 AliMillePede2.h:182
 AliMillePede2.h:183
 AliMillePede2.h:184
 AliMillePede2.h:185
 AliMillePede2.h:186
 AliMillePede2.h:187
 AliMillePede2.h:188
 AliMillePede2.h:189
 AliMillePede2.h:190
 AliMillePede2.h:191
 AliMillePede2.h:192
 AliMillePede2.h:193
 AliMillePede2.h:194
 AliMillePede2.h:195
 AliMillePede2.h:196
 AliMillePede2.h:197
 AliMillePede2.h:198
 AliMillePede2.h:199
 AliMillePede2.h:200
 AliMillePede2.h:201
 AliMillePede2.h:202
 AliMillePede2.h:203
 AliMillePede2.h:204
 AliMillePede2.h:205
 AliMillePede2.h:206
 AliMillePede2.h:207
 AliMillePede2.h:208
 AliMillePede2.h:209
 AliMillePede2.h:210
 AliMillePede2.h:211
 AliMillePede2.h:212
 AliMillePede2.h:213
 AliMillePede2.h:214
 AliMillePede2.h:215
 AliMillePede2.h:216
 AliMillePede2.h:217
 AliMillePede2.h:218
 AliMillePede2.h:219
 AliMillePede2.h:220
 AliMillePede2.h:221
 AliMillePede2.h:222
 AliMillePede2.h:223
 AliMillePede2.h:224
 AliMillePede2.h:225
 AliMillePede2.h:226
 AliMillePede2.h:227
 AliMillePede2.h:228
 AliMillePede2.h:229
 AliMillePede2.h:230
 AliMillePede2.h:231
 AliMillePede2.h:232
 AliMillePede2.h:233
 AliMillePede2.h:234
 AliMillePede2.h:235
 AliMillePede2.h:236
 AliMillePede2.h:237
 AliMillePede2.h:238
 AliMillePede2.h:239
 AliMillePede2.h:240
 AliMillePede2.h:241
 AliMillePede2.h:242
 AliMillePede2.h:243
 AliMillePede2.h:244
 AliMillePede2.h:245
 AliMillePede2.h:246
 AliMillePede2.h:247
 AliMillePede2.h:248
 AliMillePede2.h:249
 AliMillePede2.h:250
 AliMillePede2.h:251
 AliMillePede2.h:252
 AliMillePede2.h:253
 AliMillePede2.h:254
 AliMillePede2.h:255
 AliMillePede2.h:256
 AliMillePede2.h:257
 AliMillePede2.h:258
 AliMillePede2.h:259
 AliMillePede2.h:260
 AliMillePede2.h:261
 AliMillePede2.h:262
 AliMillePede2.h:263
 AliMillePede2.h:264
 AliMillePede2.h:265
 AliMillePede2.h:266
 AliMillePede2.h:267
 AliMillePede2.h:268
 AliMillePede2.h:269
 AliMillePede2.h:270
 AliMillePede2.h:271
 AliMillePede2.h:272
 AliMillePede2.h:273
 AliMillePede2.h:274
 AliMillePede2.h:275
 AliMillePede2.h:276
 AliMillePede2.h:277
 AliMillePede2.h:278
 AliMillePede2.h:279
 AliMillePede2.h:280
 AliMillePede2.h:281
 AliMillePede2.h:282
 AliMillePede2.h:283
 AliMillePede2.h:284
 AliMillePede2.h:285
 AliMillePede2.h:286
 AliMillePede2.h:287
 AliMillePede2.h:288
 AliMillePede2.h:289
 AliMillePede2.h:290
 AliMillePede2.h:291
 AliMillePede2.h:292