]> git.uio.no Git - u/mrichter/AliRoot.git/blame - STEER/AliMillePede2.h
Changing handling of the errors according to discussion following
[u/mrichter/AliRoot.git] / STEER / AliMillePede2.h
CommitLineData
8a9ab0eb 1#ifndef ALIMILLEPEDE2_H
2#define ALIMILLEPEDE2_H
3
de34b538 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
8a9ab0eb 13#include <TObject.h>
14#include <TTree.h>
8a9ab0eb 15#include "AliMinResSolve.h"
de34b538 16#include "AliMillePedeRecord.h"
17class TFile;
18class AliMatrixSq;
19class AliSymMatrix;
20class AliRectMatrix;
21class AliMatrixSparse;
8a9ab0eb 22class AliLog;
23class TStopwatch;
24
8a9ab0eb 25
8a9ab0eb 26
27class AliMillePede2: public TObject
28{
29 public:
30 //
de34b538 31 enum {gkFailed,gkInvert,gkNoInversion}; // used global matrix solution methods
8a9ab0eb 32 //
33 AliMillePede2();
34 AliMillePede2(const AliMillePede2& src);
35 virtual ~AliMillePede2();
36 AliMillePede2& operator=(const AliMillePede2& ) {printf("Dummy\n"); return *this;}
37 //
38 Int_t InitMille(int nGlo, int nLoc, Int_t lNStdDev=-1,double lResCut=-1., double lResCutInit=-1.);
39 //
40 Int_t GetNGloPar() const {return fNGloPar;}
41 Int_t GetNLocPar() const {return fNLocPar;}
42 Long_t GetNLocalEquations() const {return fNLocEquations;}
43 Int_t GetCurrentIteration() const {return fIter;}
44 Int_t GetNMaxIterations() const {return fMaxIter;}
45 Int_t GetNStdDev() const {return fNStdDev;}
46 Int_t GetNGlobalConstraints() const {return fNGloConstraints;}
de34b538 47 Int_t GetNLagrangeConstraints() const {return fNLagrangeConstraints;}
8a9ab0eb 48 Long_t GetNLocalFits() const {return fNLocFits;}
49 Long_t GetNLocalFitsRejected() const {return fNLocFitsRejected;}
50 Int_t GetNGlobalsFixed() const {return fNGloFix;}
51 Int_t GetGlobalSolveStatus() const {return fGloSolveStatus;}
52 Float_t GetChi2CutFactor() const {return fChi2CutFactor;}
53 Float_t GetChi2CutRef() const {return fChi2CutRef;}
54 Float_t GetResCurInit() const {return fResCutInit;}
55 Float_t GetResCut() const {return fResCut;}
56 Int_t GetMinPntValid() const {return fMinPntValid;}
57 Int_t GetProcessedPoints(Int_t i) const {return fProcPnt[i];}
58 Int_t* GetProcessedPoints() const {return fProcPnt;}
de34b538 59 Int_t GetParamGrID(Int_t i) const {return fParamGrID[i];}
8a9ab0eb 60 //
61 AliMatrixSq* GetGlobalMatrix() const {return fMatCGlo;}
de34b538 62 AliSymMatrix* GetLocalMatrix() const {return fMatCLoc;}
8a9ab0eb 63 Double_t* GetGlobals() const {return fVecBGlo;}
64 Double_t* GetDeltaPars() const {return fDeltaPar;}
65 Double_t* GetInitPars() const {return fInitPar;}
66 Double_t* GetSigmaPars() const {return fSigmaPar;}
de34b538 67 Bool_t* GetIsLinear() const {return fIsLinear;}
68 Double_t GetFinalParam(int i) const {return fDeltaPar[i]+fInitPar[i];}
69 Double_t GetFinalError(int i) const {return GetParError(i);}
70 //
8a9ab0eb 71 Double_t GetGlobal(Int_t i) const {return fVecBGlo[i];}
72 Double_t GetInitPar(Int_t i) const {return fInitPar[i];}
73 Double_t GetSigmaPar(Int_t i) const {return fSigmaPar[i];}
74 Bool_t GetIsLinear(Int_t i) const {return fIsLinear[i];}
75 static Bool_t IsGlobalMatSparse() {return fgIsMatGloSparse;}
de34b538 76 static Bool_t IsWeightSigma() {return fgWeightSigma;}
8a9ab0eb 77 //
de34b538 78 void SetParamGrID(Int_t grID,Int_t i) {fParamGrID[i] = grID; if(fNGroupsSet<grID)fNGroupsSet=grID;}
8a9ab0eb 79 void SetNGloPar(Int_t n) {fNGloPar = n;}
80 void SetNLocPar(Int_t n) {fNLocPar = n;}
81 void SetNMaxIterations(Int_t n=10) {fMaxIter = n;}
82 void SetNStdDev(Int_t n) {fNStdDev = n;}
83 void SetChi2CutFactor(Float_t v) {fChi2CutFactor = v;}
84 void SetChi2CutRef(Float_t v) {fChi2CutRef = v;}
85 void SetResCurInit(Float_t v) {fResCutInit = v;}
86 void SetResCut(Float_t v) {fResCut = v;}
de34b538 87 void SetMinPntValid(Int_t n) {fMinPntValid = n>0 ? n:1;}
8a9ab0eb 88 static void SetGlobalMatSparse(Bool_t v=kTRUE) {fgIsMatGloSparse = v;}
de34b538 89 static void SetWeightSigma(Bool_t v=kTRUE) {fgWeightSigma = v;}
8a9ab0eb 90 //
91 void SetInitPars(const Double_t* par) {memcpy(fInitPar,par,fNGloPar*sizeof(Double_t));}
92 void SetSigmaPars(const Double_t* par) {memcpy(fSigmaPar,par,fNGloPar*sizeof(Double_t));}
93 void SetInitPar(Int_t i,Double_t par) {fInitPar[i] = par;;}
94 void SetSigmaPar(Int_t i,Double_t par) {fSigmaPar[i] = par;}
95 //
96 Int_t GlobalFit(Double_t *par=0, Double_t *error=0, Double_t *pull=0);
97 Int_t GlobalFitIteration();
98 Int_t SolveGlobalMatEq();
de34b538 99 static void SetInvChol(Bool_t v=kTRUE) {fgInvChol = v;}
8a9ab0eb 100 static void SetMinResPrecondType(Int_t tp=0) {fgMinResCondType = tp;}
101 static void SetMinResTol(Double_t val=1e-12) {fgMinResTol = val;}
102 static void SetMinResMaxIter(Int_t val=2000) {fgMinResMaxIter = val;}
103 static void SetIterSolverType(Int_t val=AliMinResSolve::kSolMinRes) {fgIterSol = val;}
104 static void SetNKrylovV(Int_t val=60) {fgNKrylovV = val;}
105 //
de34b538 106 static Bool_t GetInvChol() {return fgInvChol;}
8a9ab0eb 107 static Int_t GetMinResPrecondType() {return fgMinResCondType;}
108 static Double_t GetMinResTol() {return fgMinResTol;}
109 static Int_t GetMinResMaxIter() {return fgMinResMaxIter;}
110 static Int_t GetIterSolverType() {return fgIterSol;}
111 static Int_t GetNKrylovV() {return fgNKrylovV;}
112 //
113 Double_t GetParError(int iPar) const;
114 Int_t PrintGlobalParameters() const;
115 //
116 //
117 Int_t SetIterations(double lChi2CutFac);
118
119 //
120 // constraints
de34b538 121 void SetGlobalConstraint(double *dergb, double val, double sigma=0);
122 void SetGlobalConstraint(const int *indgb,double *dergb, int ngb, double val, double sigma=0);
8a9ab0eb 123 //
124 // processing of the local measurement
125 void SetLocalEquation(double *dergb, double *derlc, double lMeas, double lSigma);
126 void SetLocalEquation(int *indgb, double *dergb, int ngb, int *indlc,
127 double *derlc,int nlc,double lMeas,double lSigma);
128 //
129 // manipilation with processed data and costraints records and its buffer
130 void SetDataRecFName(const char* flname) {fDataRecFName = flname;}
131 const Char_t* GetDataRecFName() const {return fDataRecFName.Data();}
132 void SetConsRecFName(const char* flname) {fConstrRecFName = flname;}
133 const Char_t* GetConsRecFName() const {return fConstrRecFName.Data();}
134
135 Bool_t InitDataRecStorage(Bool_t read=kFALSE);
136 Bool_t InitConsRecStorage(Bool_t read=kFALSE);
137 Bool_t ImposeDataRecFile(const char* fname);
138 Bool_t ImposeConsRecFile(const char* fname);
139 void CloseDataRecStorage();
140 void CloseConsRecStorage();
141 void ReadRecordData(Long_t recID);
142 void ReadRecordConstraint(Long_t recID);
143 Bool_t ReadNextRecordData();
144 Bool_t ReadNextRecordConstraint();
145 void SaveRecordData();
146 void SaveRecordConstraint();
8a9ab0eb 147 AliMillePedeRecord* GetRecord() const {return fRecord;}
148 //
149 Float_t Chi2DoFLim(int nSig, int nDoF) const;
150 //
151 // aliases for compatibility with millipede1
152 void SetParSigma(Int_t i,Double_t par) {SetSigmaPar(i,par);}
153 void SetGlobalParameters(Double_t *par) {SetInitPars(par);}
154 //
155 protected:
156 //
157 Int_t LocalFit(double *localParams=0);
158 //
159 protected:
160 //
161 Int_t fNLocPar; // number of local parameters
162 Int_t fNGloPar; // number of global parameters
163 Int_t fNGloSize; // final size of the global matrix (NGloPar+NConstraints)
164 //
165 Long_t fNLocEquations; // Number of local equations
166 Int_t fIter; // Current iteration
167 Int_t fMaxIter; // Maximum number of iterations
168 Int_t fNStdDev; // Number of standard deviations for chi2 cut
169 Int_t fNGloConstraints; // Number of constraint equations
de34b538 170 Int_t fNLagrangeConstraints; // Number of constraint equations requiring Lagrange multiplier
8a9ab0eb 171 Long_t fNLocFits; // Number of local fits
172 Long_t fNLocFitsRejected; // Number of local fits rejected
173 Int_t fNGloFix; // Number of globals fixed by user
174 Int_t fGloSolveStatus; // Status of global solver at current step
175 //
176 Float_t fChi2CutFactor; // Cut factor for chi2 cut to accept local fit
177 Float_t fChi2CutRef; // Reference cut for chi2 cut to accept local fit
178 Float_t fResCutInit; // Cut in residual for first iterartion
179 Float_t fResCut; // Cut in residual for other iterartiona
180 Int_t fMinPntValid; // min number of points for global to vary
181 //
de34b538 182 Int_t fNGroupsSet; // number of groups set
183 Int_t *fParamGrID; // group id for the every parameter
8a9ab0eb 184 Int_t *fProcPnt; // N of processed points per global variable
185 Double_t *fVecBLoc; // Vector B local (parameters)
186 Double_t *fDiagCGlo; // Initial diagonal elements of C global matrix
187 Double_t *fVecBGlo; // Vector B global (parameters)
188 //
189 Double_t *fInitPar; // Initial global parameters
190 Double_t *fDeltaPar; // Variation of global parameters
191 Double_t *fSigmaPar; // Sigma of allowed variation of global parameter
192 //
193 Bool_t *fIsLinear; // Flag for linear parameters
194 Bool_t *fConstrUsed; // Flag for used constraints
195 //
196 Int_t *fGlo2CGlo; // global ID to compressed ID buffer
197 Int_t *fCGlo2Glo; // compressed ID to global ID buffer
198 //
199 // Matrices
200 AliSymMatrix *fMatCLoc; // Matrix C local
201 AliMatrixSq *fMatCGlo; // Matrix C global
de34b538 202 AliRectMatrix *fMatCGloLoc; // Rectangular matrix C g*l
203 Int_t *fFillIndex; // auxilary index array for fast matrix fill
204 Double_t *fFillValue; // auxilary value array for fast matrix fill
8a9ab0eb 205 //
206 // processed data record bufferization
207 TString fDataRecFName; // Name of File for data records
208 AliMillePedeRecord *fRecord; // Buffer of measurements records
209 TFile *fDataRecFile; // File of processed measurements records
210 TTree *fTreeData; // Tree of processed measurements records
211 //
212 TString fConstrRecFName; // Name of File for constraints records
213 TTree *fTreeConstr; // Tree of constraint records
214 TFile *fConsRecFile; // File of processed constraints records
215 Long_t fCurrRecDataID; // ID of the current data record
216 Long_t fCurrRecConstrID; // ID of the current constraint record
de34b538 217 Bool_t fLocFitAdd; // Add contribution of carrent track (and not eliminate it)
8a9ab0eb 218 //
de34b538 219 static Bool_t fgInvChol; // Invert global matrix in Cholesky solver
220 static Bool_t fgWeightSigma; // weight parameter constraint by statistics
8a9ab0eb 221 static Bool_t fgIsMatGloSparse; // Type of the global matrix (sparse ...)
222 static Int_t fgMinResCondType; // Type of the preconditioner for MinRes method
223 static Double_t fgMinResTol; // Tolerance for MinRes solution
224 static Int_t fgMinResMaxIter; // Max number of iterations for the MinRes method
225 static Int_t fgIterSol; // type of iterative solution: MinRes or FGMRES
226 static Int_t fgNKrylovV; // size of Krylov vectors buffer in FGMRES
227 //
228 ClassDef(AliMillePede2,0)
229};
230
231//_____________________________________________________________________________________________
232inline void AliMillePede2::ReadRecordData(Long_t recID) {fTreeData->GetEntry(recID); fCurrRecDataID=recID;}
233
234//_____________________________________________________________________________________________
235inline void AliMillePede2::ReadRecordConstraint(Long_t recID) {fTreeConstr->GetEntry(recID); fCurrRecConstrID=recID;}
236
237//_____________________________________________________________________________________________
238inline void AliMillePede2::SaveRecordData() {fTreeData->Fill(); fRecord->Reset(); fCurrRecDataID++;}
239
240//_____________________________________________________________________________________________
241inline void AliMillePede2::SaveRecordConstraint() {fTreeConstr->Fill(); fRecord->Reset();fCurrRecConstrID++;}
242
243#endif