CMS 3D CMS Logo

MuonResidualsFitter.h
Go to the documentation of this file.
1 #ifndef Alignment_MuonAlignmentAlgorithms_MuonResidualsFitter_H
2 #define Alignment_MuonAlignmentAlgorithms_MuonResidualsFitter_H
3 
12 #ifndef STANDALONE_FITTER
16 #endif
17 
18 #include "TMinuit.h"
19 #include "TH1F.h"
20 #include "TProfile.h"
21 #include "TF1.h"
22 #include "TMath.h"
23 #include "TRandom3.h"
24 #include "TMatrixDSym.h"
25 
26 #include <cstdio>
27 #include <iostream>
28 #include <string>
29 #include <sstream>
30 #include <map>
31 
32 // some mock classes for the case if we want to compile fitters with pure root outside of CMSSW
33 #ifdef STANDALONE_FITTER
34 
35 #include <cassert>
36 
37 class Alignable
38 {
39 public:
40  struct Surface
41  {
42  double width() {return 300.;}
43  double length() {return 300.;}
44  };
45  Surface s;
46  Surface & surface() {return s;}
47 };
48 
49 class TFileDirectory
50 {
51 public:
52  template<typename T, typename A1, typename A2, typename A3, typename A4, typename A5>
53  T * make(const A1 &a1, const A2 &a2, const A3 &a3, const A4 &a4, const A5 &a5 ) const { T * t = new T(a1, a2, a3, a4, a5); return t; }
54  template<typename T, typename A1, typename A2, typename A3, typename A4, typename A5, typename A6, typename A7>
55  T * make(const A1 &a1, const A2 &a2, const A3 &a3, const A4 &a4, const A5 &a5, const A6 &a6, const A7 &a7) const { T * t = new T(a1, a2, a3, a4, a5, a6, a7); return t; }
56  template<typename T, typename A1, typename A2, typename A3, typename A4, typename A5, typename A6, typename A7, typename A8>
57  T * make(const A1 &a1, const A2 &a2, const A3 &a3, const A4 &a4, const A5 &a5, const A6 &a6, const A7 &a7, const A8 &a8) const { T * t = new T(a1, a2, a3, a4, a5, a6, a7, a8); return t; }
58 };
59 
60 #include <exception>
61 namespace cms{
62 class Exception : public std::exception
63 {
64 public:
65  Exception(const std::string & s): name(s) {}
66  ~Exception () throw () {}
68  template <class T> Exception& operator<<(const T& what)
69  {
70  std::cout<<name<<" exception: "<<what<<std::endl;
71  return *this;
72  }
73 };
74 }// namespace cms
75 
76 #endif // STANDALONE_FITTER
77 
78 
79 
81 {
82 public:
83  enum {
88  kPureGaussian2D
89  };
90 
91  enum {
98  kAngleBfieldFitter
99  };
100 
101  enum {
108  k0100
109  };
110 
112  {
113  Bool_t is_plus; // for +- endcap
114  Bool_t is_dt; // DT or CSC
115  UChar_t station; // 8bit uint
116  Char_t ring_wheel;
117  UChar_t sector;
118  Float_t res_x;
119  Float_t res_y;
120  Float_t res_slope_x;
121  Float_t res_slope_y;
122  Float_t pos_x;
123  Float_t pos_y;
124  Float_t angle_x;
125  Float_t angle_y;
126  Float_t pz;
127  Float_t pt;
128  Char_t q;
129  Bool_t select;
130  };
131 
133  virtual ~MuonResidualsFitter();
134 
135  virtual int type() const = 0;
136  virtual int npar() = 0;
137  virtual int ndata() = 0;
138 
139  int useRes(int pattern = -1) { if (pattern>=0) m_useResiduals = pattern; return m_useResiduals; }
140  int residualsModel() const { return m_residualsModel; }
141  long numResiduals() const { return m_residuals.size(); }
142 
143  void fix(int parNum, bool dofix=true);
144  bool fixed(int parNum);
145 
146  int nfixed() { return std::count(m_fixed.begin(), m_fixed.end(), true); }
147 
148  void setPrintLevel(int printLevel) { m_printLevel = printLevel; }
149  void setStrategy(int strategy) { m_strategy = strategy; }
150 
151  void setInitialValue(int parNum, double value) { m_parNum2InitValue[parNum] = value; }
152 
153  // an array of the actual residual and associated baggage (qoverpt, trackangle, trackposition)
154  // arrays passed to fill() are "owned" by MuonResidualsFitter: MuonResidualsFitter will delete them, don't do it yourself!
155  void fill(double *residual);
156 
157  // this block of results is only valid if fit() returns true
158  // also gamma is only valid if the model is kPowerLawTails or kROOTVoigt
159  virtual bool fit(Alignable *ali) = 0;
160 
161  double value(int parNum) { assert(0 <= parNum && parNum < npar()); return m_value[parNum]; }
162  double errorerror(int parNum) { assert(0 <= parNum && parNum < npar()); return m_error[parNum]; }
163 
164  // parNum corresponds to the parameter number defined by enums in specific fitters
165  // parIdx is a continuous index in a set of parameters
166  int parNum2parIdx(int parNum) { return m_parNum2parIdx[parNum];}
167 
168  TMatrixDSym covarianceMatrix() {return m_cov;}
169  double covarianceElement(int parNum1, int parNum2);
170  TMatrixDSym correlationMatrix();
171 
172  double loglikelihood() { return m_loglikelihood; }
173 
174  long numsegments();
175 
176  virtual double sumofweights() = 0;
177 
178  // demonstration plots; return reduced chi**2
179  virtual double plot(std::string name, TFileDirectory *dir, Alignable *ali) = 0;
180 
181  void plotsimple(std::string name, TFileDirectory *dir, int which, double multiplier);
182  void plotweighted(std::string name, TFileDirectory *dir, int which, int whichredchi2, double multiplier);
183 
184 #ifdef STANDALONE_FITTER
185  Alignable m_ali;
186  TFileDirectory m_dir;
187  bool fit() { return fit(&m_ali); }
188  virtual double plot(std::string &name) { return plot(name, &m_dir, &m_ali); }
189  void plotsimple(std::string &name, int which, double multiplier) { plotsimple(name, &m_dir, which, multiplier); }
190  void plotweighted(std::string &name, int which, int whichredchi2, double multiplier) { plotweighted(name, &m_dir, which, whichredchi2, multiplier); }
191 #endif
192 
193  // I/O of temporary files for collect mode
194  void write(FILE *file, int which=0);
195  void read(FILE *file, int which=0);
196 
197  // these are for the FCN to access what it needs to
198  std::vector<double*>::const_iterator residuals_begin() const { return m_residuals.begin(); }
199  std::vector<double*>::const_iterator residuals_end() const { return m_residuals.end(); }
200 
201  void computeHistogramRangeAndBinning(int which, int &nbins, double &a, double &b);
202  void histogramChi2GaussianFit(int which, double &fit_mean, double &fit_sigma);
203  void selectPeakResiduals_simple(double nsigma, int nvar, int *vars);
204  void selectPeakResiduals(double nsigma, int nvar, int *vars);
205 
206  // void fiducialCuts(double xMin = -1000, double xMax = 1000, double yMin = -1000, double yMax = 1000, bool fidcut1=true); // "No fiducial cut"
207  // void fiducialCuts(double xMin = -80.0, double xMax = 80.0, double yMin = -80.0, double yMax = 80.0, bool fidcut1=true); // "old" fiducial cut
208  void fiducialCuts(double xMin = -80.0, double xMax = 80.0, double yMin = -80.0, double yMax = 80.0, bool fidcut1=false); // "new" fiducial cut
209 
210  virtual void correctBField() = 0;
211  virtual void correctBField(int idx_momentum, int idx_q);
212 
213  std::vector<bool> & selectedResidualsFlags() {return m_residuals_ok;}
214 
215  void eraseNotSelectedResiduals();
216 
217 protected:
218  void initialize_table();
219  bool dofit(void (*fcn)(int&,double*,double&,double*,int), std::vector<int> &parNum, std::vector<std::string> &parName, std::vector<double> &start, std::vector<double> &step, std::vector<double> &low, std::vector<double> &high);
220  virtual void inform(TMinuit *tMinuit) = 0;
221 
226  std::vector<bool> m_fixed;
227  std::map<int, double> m_parNum2InitValue;
228  int m_printLevel, m_strategy;
229 
230  std::vector<double*> m_residuals;
231  std::vector<bool> m_residuals_ok;
232 
233  std::vector<double> m_value;
234  std::vector<double> m_error;
235  TMatrixDSym m_cov;
237 
238  std::map<int,int> m_parNum2parIdx;
239 
240  // center and radii of ellipsoid for peak selection
241  // NOTE: 10 is a hardcoded maximum of ellipsoid variables!!!
242  // but I can't imagine it ever becoming larger
243  double m_center[20];
244  double m_radii[20];
245 };
246 
247 // A ROOT-sponsored hack to get information into the fit function
248 class MuonResidualsFitterFitInfo: public TObject
249 {
250 public:
251  MuonResidualsFitterFitInfo(MuonResidualsFitter *fitter) : m_fitter(fitter) {}
252  MuonResidualsFitter *fitter() { return m_fitter; }
253 private:
255 #ifdef STANDALONE_FITTER
256  ClassDef(MuonResidualsFitterFitInfo,1);
257 #endif
258 };
259 #ifdef STANDALONE_FITTER
261 #endif
262 
263 // fit functions (these can't be put in the class; "MuonResidualsFitter_" prefix avoids namespace clashes)
264 double MuonResidualsFitter_integrate_pureGaussian(double low, double high, double center, double sigma);
265 double MuonResidualsFitter_logPureGaussian(double residual, double center, double sigma);
266 Double_t MuonResidualsFitter_pureGaussian_TF1(Double_t *xvec, Double_t *par);
267 double MuonResidualsFitter_logPureGaussian2D(double x, double y, double x0, double y0, double sx, double sy, double r);
268 double MuonResidualsFitter_compute_log_convolution(double toversigma, double gammaoversigma, double max=1000., double step=0.001, double power=4.);
269 double MuonResidualsFitter_logPowerLawTails(double residual, double center, double sigma, double gamma);
270 Double_t MuonResidualsFitter_powerLawTails_TF1(Double_t *xvec, Double_t *par);
271 double MuonResidualsFitter_logROOTVoigt(double residual, double center, double sigma, double gamma);
272 Double_t MuonResidualsFitter_ROOTVoigt_TF1(Double_t *xvec, Double_t *par);
273 double MuonResidualsFitter_logGaussPowerTails(double residual, double center, double sigma);
274 Double_t MuonResidualsFitter_GaussPowerTails_TF1(Double_t *xvec, Double_t *par);
275 void MuonResidualsPositionFitter_FCN(int &npar, double *gin, double &fval, double *par, int iflag);
276 void MuonResidualsAngleFitter_FCN(int &npar, double *gin, double &fval, double *par, int iflag);
277 
278 #endif // Alignment_MuonAlignmentAlgorithms_MuonResidualsFitter_H
Definition: start.py:1
type
Definition: HCALResponse.h:21
double MuonResidualsFitter_logGaussPowerTails(double residual, double center, double sigma)
double MuonResidualsFitter_integrate_pureGaussian(double low, double high, double center, double sigma)
residualsModel
Definition: align_cfg.py:34
Double_t MuonResidualsFitter_powerLawTails_TF1(Double_t *xvec, Double_t *par)
double MuonResidualsFitter_logPureGaussian2D(double x, double y, double x0, double y0, double sx, double sy, double r)
int parNum2parIdx(int parNum)
std::vector< double > m_error
S make(const edm::ParameterSet &cfg)
void MuonResidualsPositionFitter_FCN(int &npar, double *gin, double &fval, double *par, int iflag)
MuonResidualsFitter * fitter()
double value(int parNum)
std::ostream & operator<<(std::ostream &out, const ALILine &li)
Definition: ALILine.cc:188
TMatrixDSym covarianceMatrix()
double MuonResidualsFitter_logROOTVoigt(double residual, double center, double sigma, double gamma)
double MuonResidualsFitter_compute_log_convolution(double toversigma, double gammaoversigma, double max=1000., double step=0.001, double power=4.)
std::vector< double > m_value
tuple weightAlignment
Definition: align_cfg.py:30
std::vector< double * > m_residuals
ClassImp(AliDaqEventHeader) ClassImp(AliDaqPosition2D) ClassImp(AliDaqPositionCOPS) ClassImp(AliDaqTilt) ClassImp(AliDaqDistance) ClassImp(AliDaqTemperature) ClassImp(CocoaDaqRootEvent) CocoaDaqRootEvent
double MuonResidualsFitter_logPowerLawTails(double residual, double center, double sigma, double gamma)
std::vector< bool > m_residuals_ok
Double_t MuonResidualsFitter_GaussPowerTails_TF1(Double_t *xvec, Double_t *par)
std::map< int, int > m_parNum2parIdx
Definition: value.py:1
const AlignableSurface & surface() const
Return the Surface (global position and orientation) of the object.
Definition: Alignable.h:135
Double_t MuonResidualsFitter_pureGaussian_TF1(Double_t *xvec, Double_t *par)
Namespace of DDCMS conversion namespace.
void setStrategy(int strategy)
MuonResidualsFitter * m_fitter
double errorerror(int parNum)
std::map< int, double > m_parNum2InitValue
void fcn(int &, double *, double &, double *, int)
double b
Definition: hdecay.h:120
std::vector< double * >::const_iterator residuals_end() const
void setInitialValue(int parNum, double value)
int useRes(int pattern=-1)
std::vector< bool > & selectedResidualsFlags()
double a
Definition: hdecay.h:121
void MuonResidualsAngleFitter_FCN(int &npar, double *gin, double &fval, double *par, int iflag)
def which(cmd)
Definition: eostools.py:336
step
std::vector< bool > m_fixed
dbl *** dir
Definition: mlp_gen.cc:35
void setPrintLevel(int printLevel)
def write(self, setup)
Double_t MuonResidualsFitter_ROOTVoigt_TF1(Double_t *xvec, Double_t *par)
std::vector< double * >::const_iterator residuals_begin() const
long double T
MuonResidualsFitterFitInfo(MuonResidualsFitter *fitter)
double MuonResidualsFitter_logPureGaussian(double residual, double center, double sigma)