]> git.uio.no Git - u/mrichter/AliRoot.git/blob - STEER/STEER/AliMillePede2.h
a91c531595196e9b1e5f3709bbe337575bbb4f6c
[u/mrichter/AliRoot.git] / STEER / STEER / AliMillePede2.h
1 #ifndef ALIMILLEPEDE2_H
2 #define ALIMILLEPEDE2_H
3
4 /**********************************************************************************************/
5 /* General class for alignment with large number of degrees of freedom                        */
6 /* Based on the original milliped2 by Volker Blobel                                           */
7 /* http://www.desy.de/~blobel/mptalks.html                                                    */
8 /*                                                                                            */ 
9 /* Author: ruben.shahoyan@cern.ch                                                             */
10 /*                                                                                            */ 
11 /**********************************************************************************************/
12
13 #include <TObject.h>
14 #include <TString.h>
15 #include <TTree.h>
16 #include "AliMinResSolve.h"
17 #include "AliMillePedeRecord.h"
18 class TFile;
19 class AliMatrixSq;
20 class AliSymMatrix;
21 class AliRectMatrix;
22 class AliMatrixSparse;
23 class AliLog;
24 class TStopwatch;
25 class TArrayL;
26 class TArrayF;
27
28
29 class AliMillePede2: public TObject
30 {
31  public:
32   //
33   enum {kFailed,kInvert,kNoInversion};    // used global matrix solution methods
34   enum {kFixParID=-1};                    // dummy id for fixed param
35   //
36   AliMillePede2();
37   AliMillePede2(const AliMillePede2& src);
38   virtual ~AliMillePede2();
39   AliMillePede2& operator=(const AliMillePede2& )             {printf("Dummy\n"); return *this;}
40   //
41   Int_t                InitMille(int nGlo, int nLoc, Int_t lNStdDev=-1,double lResCut=-1., double lResCutInit=-1., const Int_t* regroup=0);
42   //
43   Int_t                GetNGloPar()                     const {return fNGloPar;}
44   Int_t                GetNGloParIni()                  const {return fNGloParIni;}
45   const Int_t*         GetRegrouping()                  const {return fkReGroup;}
46   Int_t                GetNLocPar()                     const {return fNLocPar;}
47   Long_t               GetNLocalEquations()             const {return fNLocEquations;}
48   Int_t                GetCurrentIteration()            const {return fIter;}
49   Int_t                GetNMaxIterations()              const {return fMaxIter;}
50   Int_t                GetNStdDev()                     const {return fNStdDev;} 
51   Int_t                GetNGlobalConstraints()          const {return fNGloConstraints;}
52   Int_t                GetNLagrangeConstraints()        const {return fNLagrangeConstraints;}
53   Long_t               GetNLocalFits()                  const {return fNLocFits;}
54   Long_t               GetNLocalFitsRejected()          const {return fNLocFitsRejected;}
55   Int_t                GetNGlobalsFixed()               const {return fNGloFix;}
56   Int_t                GetGlobalSolveStatus()           const {return fGloSolveStatus;}
57   Float_t              GetChi2CutFactor()               const {return fChi2CutFactor;}
58   Float_t              GetChi2CutRef()                  const {return fChi2CutRef;}
59   Float_t              GetResCurInit()                  const {return fResCutInit;}
60   Float_t              GetResCut()                      const {return fResCut;}
61   Int_t                GetMinPntValid()                 const {return fMinPntValid;}
62   Int_t                GetRGId(Int_t i)                 const {return fkReGroup ? (fkReGroup[i]<0 ? -1:fkReGroup[i]) : i;}       
63   Int_t                GetProcessedPoints(Int_t i)      const {int ir=GetRGId(i); return ir<=0 ? 0:fProcPnt[ir];}
64   Int_t*               GetProcessedPoints()             const {return fProcPnt;}
65   Int_t                GetParamGrID(Int_t i)            const {int ir=GetRGId(i); return ir<=0 ? 0:fParamGrID[ir];}
66   //
67   AliMatrixSq*         GetGlobalMatrix()                const {return fMatCGlo;}
68   AliSymMatrix*        GetLocalMatrix()                 const {return fMatCLoc;}
69   Double_t*            GetGlobals()                     const {return fVecBGlo;}
70   Double_t*            GetDeltaPars()                   const {return fDeltaPar;}
71   Double_t*            GetInitPars()                    const {return fInitPar;}
72   Double_t*            GetSigmaPars()                   const {return fSigmaPar;}
73   Bool_t*              GetIsLinear()                    const {return fIsLinear;}
74   Double_t             GetFinalParam(int i)             const {int ir=GetRGId(i); return ir<0 ? 0:fDeltaPar[ir]+fInitPar[ir];}
75   Double_t             GetFinalError(int i)             const {return GetParError(i);}
76   Double_t             GetPull(int i)                   const;
77   //
78   Double_t             GetGlobal(Int_t i)               const {int ir=GetRGId(i); return ir<0 ? 0:fVecBGlo[ir];}
79   Double_t             GetInitPar(Int_t i)              const {int ir=GetRGId(i); return ir<0 ? 0:fInitPar[ir];}
80   Double_t             GetSigmaPar(Int_t i)             const {int ir=GetRGId(i); return ir<0 ? 0:fSigmaPar[ir];}
81   Bool_t               GetIsLinear(Int_t i)             const {int ir=GetRGId(i); return ir<0 ? 0:fIsLinear[ir];}
82   static Bool_t        IsGlobalMatSparse()                    {return fgIsMatGloSparse;}
83   static Bool_t        IsWeightSigma()                        {return fgWeightSigma;}
84   //
85   void                 SetParamGrID(Int_t grID,Int_t i)       {int ir=GetRGId(i); if(ir<0) return; fParamGrID[ir] = grID; if(fNGroupsSet<grID)fNGroupsSet=grID;}
86   void                 SetNGloPar(Int_t n)                    {fNGloPar = n;}
87   void                 SetNLocPar(Int_t n)                    {fNLocPar = n;}
88   void                 SetNMaxIterations(Int_t n=10)          {fMaxIter = n;}
89   void                 SetNStdDev(Int_t n)                    {fNStdDev = n;}
90   void                 SetChi2CutFactor(Float_t v)            {fChi2CutFactor = v;}
91   void                 SetChi2CutRef(Float_t v)               {fChi2CutRef = v;}
92   void                 SetResCurInit(Float_t v)               {fResCutInit = v;}
93   void                 SetResCut(Float_t v)                   {fResCut = v;}
94   void                 SetMinPntValid(Int_t n)                {fMinPntValid = n>0 ? n:1;}
95   static void          SetGlobalMatSparse(Bool_t v=kTRUE)     {fgIsMatGloSparse = v;}
96   static void          SetWeightSigma(Bool_t v=kTRUE)         {fgWeightSigma = v;}
97   //
98   void                 SetInitPars(const Double_t* par);
99   void                 SetSigmaPars(const Double_t* par);
100   void                 SetInitPar(Int_t i,Double_t par);
101   void                 SetSigmaPar(Int_t i,Double_t par);
102   //
103   Int_t                GlobalFit(Double_t *par=0, Double_t *error=0, Double_t *pull=0);
104   Int_t                GlobalFitIteration();
105   Int_t                SolveGlobalMatEq();
106   static void          SetInvChol(Bool_t v=kTRUE)             {fgInvChol = v;}
107   static void          SetMinResPrecondType(Int_t tp=0)       {fgMinResCondType = tp;}
108   static void          SetMinResTol(Double_t val=1e-12)       {fgMinResTol = val;}
109   static void          SetMinResMaxIter(Int_t val=2000)       {fgMinResMaxIter = val;}
110   static void          SetIterSolverType(Int_t val=AliMinResSolve::kSolMinRes) {fgIterSol = val;}
111   static void          SetNKrylovV(Int_t val=60)              {fgNKrylovV = val;}
112   //
113   static Bool_t        GetInvChol()                           {return fgInvChol;}
114   static Int_t         GetMinResPrecondType()                 {return fgMinResCondType;}
115   static Double_t      GetMinResTol()                         {return fgMinResTol;}
116   static Int_t         GetMinResMaxIter()                     {return fgMinResMaxIter;}
117   static Int_t         GetIterSolverType()                    {return fgIterSol;}
118   static Int_t         GetNKrylovV()                          {return fgNKrylovV;}
119   //
120   Double_t             GetParError(int iPar)           const;
121   Int_t                PrintGlobalParameters()         const;
122   void                 SetRejRunList(const UInt_t *runs, Int_t nruns);
123   void                 SetAccRunList(const UInt_t *runs, Int_t nruns, const Float_t* wghList=0);
124   Bool_t               IsRecordAcceptable();
125   //
126   //
127   Int_t                SetIterations(double lChi2CutFac);
128
129   //
130   // constraints
131   void                 SetGlobalConstraint(const double *dergb, double val, double sigma=0);
132   void                 SetGlobalConstraint(const int *indgb, const double *dergb, int ngb, double val, double sigma=0);
133   //
134   // processing of the local measurement
135   void                 SetRecordRun(Int_t run);
136   void                 SetRecordWeight(double wgh);
137   void                 SetLocalEquation(double *dergb, double *derlc, double lMeas, double lSigma);
138   void                 SetLocalEquation(int *indgb, double *dergb, int ngb, int *indlc, 
139                                         double *derlc,int nlc,double lMeas,double lSigma);
140   //
141   // manipilation with processed data and costraints records and its buffer
142   void                 SetDataRecFName(const char* flname)   {fDataRecFName = flname;}
143   const Char_t*        GetDataRecFName()               const {return fDataRecFName.Data();}
144   void                 SetConsRecFName(const char* flname)   {fConstrRecFName = flname;}
145   const Char_t*        GetConsRecFName()               const {return fConstrRecFName.Data();}
146   //
147   void   SetRecDataTreeName(const char* name=0)     {fRecDataTreeName = name;   if (fRecDataTreeName.IsNull()) fRecDataTreeName = "AliMillePedeRecords_Data";}
148   void   SetRecConsTreeName(const char* name=0)     {fRecConsTreeName = name;   if (fRecConsTreeName.IsNull()) fRecConsTreeName = "AliMillePedeRecords_Consaints";}
149   void   SetRecDataBranchName(const char* name=0)   {fRecDataBranchName = name; if (fRecDataBranchName.IsNull()) fRecDataBranchName = "Record_Data";}
150   void   SetRecConsBranchName(const char* name=0)   {fRecConsBranchName = name; if (fRecConsBranchName.IsNull()) fRecConsBranchName = "Record_Consaints";}
151   const char* GetRecDataTreeName()     const {return fRecDataTreeName.Data();}
152   const char* GetRecConsTreeName()     const {return fRecConsTreeName.Data();}
153   const char* GetRecDataBranchName()   const {return fRecDataBranchName.Data();}
154   const char* GetRecConsBranchName()   const {return fRecConsBranchName.Data();}
155   //
156   Bool_t               InitDataRecStorage(Bool_t read=kFALSE);
157   Bool_t               InitConsRecStorage(Bool_t read=kFALSE);
158   Bool_t               ImposeDataRecFile(const char* fname);
159   Bool_t               ImposeConsRecFile(const char* fname);
160   void                 CloseDataRecStorage();
161   void                 CloseConsRecStorage();
162   void                 ReadRecordData(Long_t recID);
163   void                 ReadRecordConstraint(Long_t recID);
164   Bool_t               ReadNextRecordData();
165   Bool_t               ReadNextRecordConstraint();
166   void                 SaveRecordData();
167   void                 SaveRecordConstraint();
168   AliMillePedeRecord*  GetRecord()                      const {return fRecord;}
169   Long_t               GetSelFirst()                    const {return fSelFirst;}
170   Long_t               GetSelLast()                     const {return fSelLast;}
171   void                 SetSelFirst(Long_t v)                  {fSelFirst = v;}
172   void                 SetSelLast(Long_t v)                   {fSelLast = v;}
173   //
174   Float_t              Chi2DoFLim(int nSig, int nDoF)   const;
175   //
176   // aliases for compatibility with millipede1
177   void                 SetParSigma(Int_t i,Double_t par)      {SetSigmaPar(i,par);}
178   void                 SetGlobalParameters(Double_t *par)     {SetInitPars(par);}
179   void                 SetNonLinear(int index, Bool_t v=kTRUE) {int id = GetRGId(index); if (id<0) return; fIsLinear[id] = !v;}
180   //
181  protected:
182   //
183   Int_t                LocalFit(double *localParams=0);
184   Bool_t               IsZero(Double_t v,Double_t eps=1e-16)   const {return TMath::Abs(v)<eps;}                  
185   //
186  protected:
187   //
188   Int_t                 fNLocPar;                        // number of local parameters
189   Int_t                 fNGloPar;                        // number of global parameters
190   Int_t                 fNGloParIni;                     // number of global parameters before grouping
191   Int_t                 fNGloSize;                       // final size of the global matrix (NGloPar+NConstraints)
192   //
193   Long_t                fNLocEquations;                  // Number of local equations 
194   Int_t                 fIter;                           // Current iteration
195   Int_t                 fMaxIter;                        // Maximum number of iterations
196   Int_t                 fNStdDev;                        // Number of standard deviations for chi2 cut 
197   Int_t                 fNGloConstraints;                // Number of constraint equations
198   Int_t                 fNLagrangeConstraints;           // Number of constraint equations requiring Lagrange multiplier
199   Long_t                fNLocFits;                       // Number of local fits
200   Long_t                fNLocFitsRejected;               // Number of local fits rejected
201   Int_t                 fNGloFix;                        // Number of globals fixed by user
202   Int_t                 fGloSolveStatus;                 // Status of global solver at current step
203   //
204   Float_t               fChi2CutFactor;                  // Cut factor for chi2 cut to accept local fit 
205   Float_t               fChi2CutRef;                     // Reference cut for chi2 cut to accept local fit 
206   Float_t               fResCutInit;                     // Cut in residual for first iterartion
207   Float_t               fResCut;                         // Cut in residual for other iterartiona
208   Int_t                 fMinPntValid;                    // min number of points for global to vary
209   //
210   Int_t                 fNGroupsSet;                     // number of groups set
211   Int_t                *fParamGrID;                      //[fNGloPar] group id for the every parameter
212   Int_t                *fProcPnt;                        //[fNGloPar] N of processed points per global variable
213   Double_t             *fVecBLoc;                        //[fNLocPar] Vector B local (parameters) 
214   Double_t             *fDiagCGlo;                       //[fNGloPar] Initial diagonal elements of C global matrix
215   Double_t             *fVecBGlo;                        //! Vector B global (parameters)
216   //
217   Double_t             *fInitPar;                        //[fNGloPar] Initial global parameters
218   Double_t             *fDeltaPar;                       //[fNGloPar] Variation of global parameters
219   Double_t             *fSigmaPar;                       //[fNGloPar] Sigma of allowed variation of global parameter
220   //
221   Bool_t               *fIsLinear;                       //[fNGloPar] Flag for linear parameters
222   Bool_t               *fConstrUsed;                     //! Flag for used constraints
223   //
224   Int_t                *fGlo2CGlo;                       //[fNGloPar] global ID to compressed ID buffer
225   Int_t                *fCGlo2Glo;                       //[fNGloPar] compressed ID to global ID buffer
226   //
227   // Matrices
228   AliSymMatrix         *fMatCLoc;                        // Matrix C local
229   AliMatrixSq          *fMatCGlo;                        // Matrix C global
230   AliRectMatrix        *fMatCGloLoc;                     // Rectangular matrix C g*l 
231   Int_t                *fFillIndex;                      //[fNGloPar] auxilary index array for fast matrix fill
232   Double_t             *fFillValue;                      //[fNGloPar] auxilary value array for fast matrix fill
233   //
234   // processed data record bufferization   
235   TString               fRecDataTreeName;                // Name of data records tree
236   TString               fRecConsTreeName;                // Name of constraints records tree
237   TString               fRecDataBranchName;              // Name of data records branch name
238   TString               fRecConsBranchName;              // Name of constraints records branch name
239   
240   TString               fDataRecFName;                   // Name of File for data records               
241   AliMillePedeRecord   *fRecord;                         // Buffer of measurements records
242   TFile                *fDataRecFile;                    // File of processed measurements records
243   TTree                *fTreeData;                       // Tree of processed measurements records
244   Int_t                 fRecFileStatus;                  // state of the record file (0-no, 1-read, 2-rw)
245   //
246   TString               fConstrRecFName;                 // Name of File for constraints records               
247   TTree                *fTreeConstr;                     //! Tree of constraint records
248   TFile                *fConsRecFile;                    //! File of processed constraints records
249   Long_t                fCurrRecDataID;                  // ID of the current data record
250   Long_t                fCurrRecConstrID;                // ID of the current constraint record
251   Bool_t                fLocFitAdd;                      // Add contribution of carrent track (and not eliminate it)
252   Int_t                 fSelFirst;                       // event selection start
253   Int_t                 fSelLast;                        // event selection end
254   TArrayL*              fRejRunList;                     // list of runs to reject (if any)
255   TArrayL*              fAccRunList;                     // list of runs to select (if any)
256   TArrayF*              fAccRunListWgh;                  // optional weights for data of accepted runs (if any)
257   Double_t              fRunWgh;                         // run weight
258   const Int_t*          fkReGroup;                       // optional regrouping of parameters wrt ID's from the records
259   //
260   static Bool_t         fgInvChol;                       // Invert global matrix in Cholesky solver
261   static Bool_t         fgWeightSigma;                   // weight parameter constraint by statistics
262   static Bool_t         fgIsMatGloSparse;                // Type of the global matrix (sparse ...)
263   static Int_t          fgMinResCondType;                // Type of the preconditioner for MinRes method 
264   static Double_t       fgMinResTol;                     // Tolerance for MinRes solution
265   static Int_t          fgMinResMaxIter;                 // Max number of iterations for the MinRes method
266   static Int_t          fgIterSol;                       // type of iterative solution: MinRes or FGMRES
267   static Int_t          fgNKrylovV;                      // size of Krylov vectors buffer in FGMRES
268   //
269   ClassDef(AliMillePede2,1)
270 };
271
272 //_____________________________________________________________________________________________
273 inline void AliMillePede2::ReadRecordData(Long_t recID)       {fTreeData->GetEntry(recID); fCurrRecDataID=recID;}
274
275 //_____________________________________________________________________________________________
276 inline void AliMillePede2::ReadRecordConstraint(Long_t recID) {fTreeConstr->GetEntry(recID); fCurrRecConstrID=recID;}
277
278 //_____________________________________________________________________________________________
279 inline void AliMillePede2::SaveRecordData()                   {fTreeData->Fill(); fRecord->Reset(); fCurrRecDataID++;}
280
281 //_____________________________________________________________________________________________
282 inline void AliMillePede2::SaveRecordConstraint()             {fTreeConstr->Fill(); fRecord->Reset();fCurrRecConstrID++;}
283
284 #endif