ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/JetFitAnalyzer/src/jetfit.cpp
Revision: 1.12
Committed: Wed Nov 11 23:20:31 2009 UTC (15 years, 5 months ago) by dnisson
Branch: MAIN
Changes since 1.11: +2 -164 lines
Log Message:
Removed erfc and prob

File Contents

# User Rev Content
1 dnisson 1.1 /* jetfit.cpp - Package to fit multi-Gaussian distributions to histograms
2     * rewrite after accidental deletion 07-24-09
3     * Author: David Nisson
4     */
5    
6     #include <cstdlib>
7     #include <cmath>
8     #include <ctime>
9     #include <cstdio>
10     #include <cstring>
11     #include <cctype>
12     #include <iostream>
13     #include <fstream>
14     #include <sstream>
15     #include <string>
16     #include <vector>
17     #include <set>
18     #include <map>
19     #include <limits>
20     #include <exception>
21    
22     #include <TMinuit.h>
23     #include <TH1.h>
24     #include <TH2.h>
25     #include <TF2.h>
26     #include <TRandom3.h>
27     #include <Rtypes.h>
28     #include <TFormula.h>
29     #include <TSystem.h>
30 dnisson 1.12 #include <TMath.h>
31 dnisson 1.1
32     #include "UserCode/JetFitAnalyzer/interface/jetfit.h"
33    
34     #define PI 3.141592653
35    
36     #define MAX_GAUSS 3
37    
38     using namespace std;
39    
40     namespace jetfit {
41     // configurable options
42     bool ignorezero = false; // ignore zero bins when fitting
43    
44     struct histo {
45     vector< vector<double> > bins;
46     double Xlo, Xhi, Ylo, Yhi;
47     } energy;
48    
49    
50     model_def *mdef;
51    
52     model_def& curr_model_def() {
53     return *mdef;
54     }
55    
56     void set_model_def(model_def *_mdef) {
57     mdef = _mdef;
58     }
59    
60 dnisson 1.10 int model_def::get_ngauss() {
61 dnisson 1.1 return ngauss;
62     }
63    
64 dnisson 1.10 void model_def::set_ngauss(int _ngauss) {
65 dnisson 1.1 ngauss = _ngauss;
66     }
67    
68     // fit function
69 dnisson 1.10 double model_def::fit_fcn(double x, double y, double *xval) {
70     int npar_indiv = get_formula()->GetNpar();
71 dnisson 1.1 double val = 0.0;
72     for (int i = 0; i < ngauss; i++) {
73 dnisson 1.10 get_formula()->SetParameters(xval + i*npar_indiv);
74 dnisson 1.1 val += mdef->get_formula()->Eval(x, y);
75     }
76     return val;
77     }
78    
79     // TF2-compatible fit function
80     double fit_fcn_TF2(double *x, double *pval) {
81 dnisson 1.10 double val = mdef->fit_fcn(x[0], x[1], pval);
82 dnisson 1.1 return val;
83     }
84    
85 dnisson 1.9 // Integral of (model formula)^2 / chisquare sigma
86 dnisson 1.10 double model_def::formula_int(double xlo, double xhi, double ylo, double yhi,
87 dnisson 1.1 double *pval, double XbinSize, double YbinSize,
88     double *xval) {
89     double xstep = (xhi - xlo) / static_cast<double>(20);
90     double ystep = (yhi - ylo) / static_cast<double>(20);
91    
92     double fsum = 0.0;
93     double pval_old[256];
94 dnisson 1.10 int npar = get_formula()->GetNpar();
95 dnisson 1.1 if (npar > 256) {
96     cerr << "Parameter overload" << endl;
97     return 0.0;
98     }
99 dnisson 1.10 get_formula()->GetParameters(pval_old);
100     get_formula()->SetParameters(pval);
101 dnisson 1.1
102     for (int i = 0; i < 20; i++) {
103     for (int j = 0; j < 20; j++) {
104     double x = (static_cast<double>(i) + 0.5)*xstep + xlo;
105     double y = (static_cast<double>(j) + 0.5)*ystep + ylo;
106     double sig_cut = 1.0e-5;
107     double nu = XbinSize * YbinSize * fit_fcn(x, y, xval);
108 dnisson 1.10 double sigma = fabs(chisquare_error(0.0)) < sig_cut
109     ? sig_cut : chisquare_error(0.0);
110 dnisson 1.1
111 dnisson 1.10 fsum += pow(xstep * ystep * get_formula()->Eval(x, y), 2.0)
112 dnisson 1.1 / sigma / sigma;
113     }
114     }
115    
116 dnisson 1.10 get_formula()->SetParameters(pval_old);
117 dnisson 1.1 return fsum;
118     }
119    
120     // MINUIT objective function
121     void fcn(int &npar, double *grad, double &fcnval, double *xval, int iflag)
122     {
123     double chisquare = 0.0;
124     double XbinSize = (energy.Xhi - energy.Xlo)
125     / static_cast<double>(energy.bins.size());
126     double YbinSize = (energy.Yhi - energy.Ylo)
127     / static_cast<double>(energy.bins.at(0).size());
128     try {
129     // add errors in data points in histogram
130     for (int i = 0; i < energy.bins.size(); i++) {
131     for (int j = 0; j < energy.bins.at(i).size(); j++) {
132 dnisson 1.4 double x = (static_cast<double>(i) + 0.5) * (energy.Xhi - energy.Xlo)
133 dnisson 1.1 / static_cast<double>(energy.bins.size()) + energy.Xlo;
134     double y = (static_cast<double>(j) + 0.5) * (energy.Yhi - energy.Ylo)
135     / static_cast<double>(energy.bins.at(i).size()) + energy.Ylo;
136 dnisson 1.10 double nu = mdef->fit_fcn(x, y, xval) * XbinSize * YbinSize;
137 dnisson 1.1 double sig_cut = 1.0e-5;
138 dnisson 1.7 if (fabs(mdef->chisquare_error(energy.bins.at(i).at(j))) > sig_cut || !ignorezero) {
139 dnisson 1.1 chisquare += pow(energy.bins.at(i).at(j) - nu, 2.0)
140 dnisson 1.7 / pow(mdef->chisquare_error(energy.bins.at(i).at(j)), 2.0);
141 dnisson 1.3 }
142     else {
143 dnisson 1.1 chisquare += pow(energy.bins.at(i).at(j) - nu, 2.0)
144     / pow(sig_cut, 2.0);
145 dnisson 1.3 }
146 dnisson 1.1 }
147     }
148 dnisson 1.7
149 dnisson 1.1 // add errors due to Gaussians outside histogram
150 dnisson 1.6 double eps = 0.01; // accuracy set for this function
151 dnisson 1.10 for (int i = 0; i < mdef->get_ngauss(); i++) {
152 dnisson 1.3 double *pval = xval+i*(mdef->get_formula()->GetNpar());
153 dnisson 1.1 double par_x = pval[mdef->get_indiv_max_x()];
154     double par_y = pval[mdef->get_indiv_max_y()];
155     double par_sig = pval[3];
156 dnisson 1.7 double cutoff_rad = sqrt(-log(eps * 2.0 * PI * par_sig * par_sig
157     / pval[0]) * 2.0 * par_sig * par_sig);
158 dnisson 1.1
159     bool left = par_x < energy.Xlo + cutoff_rad;
160     bool right = par_x > energy.Xhi - cutoff_rad;
161     bool top = par_y > energy.Yhi - cutoff_rad;
162 dnisson 1.2 bool bottom = par_y < energy.Ylo + cutoff_rad;
163 dnisson 1.1
164     if (left) {
165     double xlo = par_x - cutoff_rad;
166     double xhi = energy.Xlo;
167     double ylo = energy.Ylo;
168     double yhi = energy.Yhi;
169 dnisson 1.10 chisquare += 2.0*mdef->formula_int(xlo, xhi, ylo, yhi, pval,
170 dnisson 1.1 XbinSize, YbinSize, xval);
171     }
172     if (right) {
173     double xlo = energy.Xhi;
174     double xhi = par_x + cutoff_rad;
175     double ylo = energy.Ylo;
176     double yhi = energy.Yhi;
177 dnisson 1.10 chisquare += 2.0*mdef->formula_int(xlo, xhi, ylo, yhi, pval,
178 dnisson 1.1 XbinSize, YbinSize, xval);
179     }
180     if (top) {
181     double xlo = left ? par_x - cutoff_rad : energy.Xlo;
182     double xhi = right ? par_x + cutoff_rad : energy.Xhi;
183     double ylo = energy.Yhi;
184     double yhi = par_y + cutoff_rad;
185 dnisson 1.10 chisquare += 2.0*mdef->formula_int(xlo, xhi, ylo, yhi, pval,
186 dnisson 1.1 XbinSize, YbinSize, xval);
187     }
188     if (bottom) {
189     double xlo = left ? par_x - cutoff_rad : energy.Xlo;
190     double xhi = right ? par_x + cutoff_rad : energy.Xhi;
191     double ylo = par_y - cutoff_rad;
192     double yhi = energy.Ylo;
193 dnisson 1.10 chisquare += 2.0*mdef->formula_int(xlo, xhi, ylo, yhi, pval,
194 dnisson 1.1 XbinSize, YbinSize, xval);
195     }
196     }
197 dnisson 1.7
198 dnisson 1.1 fcnval = chisquare;
199     }
200     catch (exception ex) {
201     cerr << "Exception in jetfit::fcn" << endl;
202     }
203     }
204    
205     bool get_ignorezero() {
206     return ignorezero;
207     }
208    
209     void set_ignorezero(bool _iz) {
210     ignorezero = _iz;
211     }
212    
213     TFormula * model_def::get_formula() {
214     return indiv_formula;
215     }
216    
217     void model_def::set_formula(TFormula *new_formula) {
218     indiv_formula = new_formula;
219     }
220    
221     void model_def::get_indiv_par(int i, string &name, double &start, double &step,
222     double &lo, double &hi) {
223     try {
224     name = indiv_par_names.at(i);
225     start = indiv_par_starts.at(i);
226     step = indiv_par_start_steps.at(i);
227     lo = indiv_par_lo.at(i);
228     hi = indiv_par_hi.at(i);
229     } catch (exception ex) {
230     cerr << "exception in jetfit::model_def::get_indiv_par" << endl;
231     }
232     }
233    
234     void model_def::set_indiv_par(int i, string name, double start, double step,
235     double lo, double hi) {
236     if (indiv_par_names.size() < i+1) {
237     indiv_par_names.resize(i+1);
238     indiv_par_starts.resize(i+1);
239     indiv_par_start_steps.resize(i+1);
240     indiv_par_lo.resize(i+1);
241     indiv_par_hi.resize(i+1);
242     }
243     try {
244     indiv_par_names.at(i) = name;
245     indiv_par_starts.at(i) = start;
246     indiv_par_start_steps.at(i) = step;
247     indiv_par_lo.at(i) = lo;
248     indiv_par_hi.at(i) = hi;
249     }
250     catch (exception ex) {
251     cerr << "exception in jetfit::model_def::set_indiv_par" << endl;
252     }
253     }
254    
255     int model_def::get_indiv_max_E() {
256     return indiv_max_E;
257     }
258    
259     void model_def::set_indiv_max_E(int _new_max_E) {
260     indiv_max_E = _new_max_E;
261     }
262    
263     int model_def::get_indiv_max_x() {
264     return indiv_max_x;
265     }
266    
267     void model_def::set_indiv_max_x(int new_max_x) {
268     indiv_max_x = new_max_x;
269     }
270    
271     int model_def::get_indiv_max_y() {
272     return indiv_max_y;
273     }
274    
275     void model_def::set_indiv_max_y(int new_max_y) {
276     indiv_max_y = new_max_y;
277     }
278    
279     unsigned model_def::get_n_special_par_sets() {
280     return special_par_starts.size();
281     }
282    
283     void model_def::get_special_par(int g, int i, double &pstart,
284     double &perr, double &plo, double &phi) {
285     try {
286     pstart = special_par_starts.at(g).at(i);
287     perr = special_par_start_steps.at(g).at(i);
288     plo = special_par_lo.at(g).at(i);
289     phi = special_par_hi.at(g).at(i);
290     } catch (exception ex) {
291     cerr << "Exception in model_def::get_special_par" << endl;
292     }
293     }
294    
295     void model_def::set_special_par(int g, int i, double pstart,
296     double perr, double plo, double phi) {
297     if (g+1 > special_par_starts.size()) {
298     special_par_starts.resize(g+1);
299     special_par_start_steps.resize(g+1);
300     special_par_lo.resize(g+1);
301     special_par_hi.resize(g+1);
302     }
303    
304     try {
305     if (i+1 > special_par_starts.at(g).size()) {
306     special_par_starts.at(g).resize(i+1);
307     special_par_start_steps.at(g).resize(i+1);
308     special_par_lo.at(g).resize(i+1);
309     special_par_hi.at(g).resize(i+1);
310     }
311    
312     special_par_starts.at(g).at(i) = pstart;
313     special_par_start_steps.at(g).at(i) = perr;
314     special_par_lo.at(g).at(i) = plo;
315     special_par_hi.at(g).at(i) = phi;
316     }
317     catch (exception ex) {
318     cerr << "exception in jetfit::model_def::set_special_par" << endl;
319     }
320     }
321    
322 dnisson 1.5 void fill_histo_with_vec(TH2 *hist, vector< vector<double> > &vec) {
323     for (int i = 0; i < vec.size(); i++) {
324     for (int j = 0; j < vec[i].size(); j++) {
325 dnisson 1.6 hist->SetBinContent(i+1, j+1, vec[i][j]);
326 dnisson 1.5 }
327     }
328     }
329    
330 dnisson 1.10 void model_def::par_init(TH2 *hist, TMinuit *gMinuit, vector<TString> &pars,
331 dnisson 1.1 double *pval, double *perr, double *plo, double *phi,
332     int npar, results r) {
333 dnisson 1.10 int npar1 = npar - get_formula()->GetNpar();
334 dnisson 1.1 bool init_spec_pars = false;
335 dnisson 1.10 if (ngauss <= get_n_special_par_sets()) {
336 dnisson 1.1 init_spec_pars = true;
337     for (int k = 0; k < ngauss; k++) {
338 dnisson 1.10 int ipar0 = k*get_formula()->GetNpar(); // index of indiv par 0
339     for (int l = 0; l < get_formula()->GetNpar(); l++) {
340 dnisson 1.1 int ipar = ipar0 + l;
341     if (ipar < npar)
342 dnisson 1.10 get_special_par(k, l, pval[ipar], perr[ipar], plo[ipar],
343 dnisson 1.1 phi[ipar]);
344     else
345     cerr << "WARNING: Attempt to set parameter out of index range!"
346     << endl;
347     }
348     }
349     }
350     else {
351     double XbinSize = (energy.Xhi - energy.Xlo)
352     / static_cast<double>(energy.bins.size());
353     double YbinSize = (energy.Xhi - energy.Xlo)
354     / static_cast<double>(energy.bins.at(0).size());
355     vector< vector<double> > sub_energy(energy.bins.size());
356     double max_sub, max_x = 0.0, max_y = 0.0;
357     double pval_other[256];
358 dnisson 1.3 max_sub = -(numeric_limits<double>::max());
359 dnisson 1.1 for (int i = 0; i < energy.bins.size(); i++) {
360     sub_energy[i].resize(energy.bins[i].size());
361     for (int j = 0; j < energy.bins[i].size(); j++) {
362     double x = static_cast<double>(i)*(energy.Xhi - energy.Xlo)
363     / static_cast<double>(energy.bins.size()) + energy.Xlo;
364     double y = static_cast<double>(j)*(energy.Yhi - energy.Ylo)
365     / static_cast<double>(energy.bins[i].size()) + energy.Ylo;
366     if (ngauss > 1) {
367 dnisson 1.6 // subtract integral of fit function
368 dnisson 1.1 try {
369     int npval_other = r.pval.at(r.pval.size()-1).size();
370     if (npval_other > 256) {
371     cerr << "Parameter overload" << endl;
372     return;
373     }
374     for (int k = 0; k < npval_other; k++) {
375     pval_other[k] = r.pval[r.pval.size()-1][k];
376     }
377     ngauss--;
378     double nu = fit_fcn(x, y, pval_other) * XbinSize * YbinSize;
379     ngauss++;
380 dnisson 1.6 sub_energy[i][j] = energy.bins[i][j] - nu;
381 dnisson 1.1 }
382     catch (exception ex) {
383     cerr << "Exception in par_init" << endl;
384     }
385     }
386     else {
387     sub_energy[i][j] = energy.bins[i][j];
388     }
389     }
390     }
391 dnisson 1.5 TH2D *hist_sub = new TH2D("hist_sub", "Subtracted histo",
392     sub_energy.size(), energy.Xlo,
393     energy.Xhi,
394     sub_energy[0].size(), energy.Ylo,
395     energy.Yhi);
396     fill_histo_with_vec(hist_sub, sub_energy);
397 dnisson 1.8 max_sub = 0.0;
398 dnisson 1.9 max_x = 0.0;
399     max_y = 0.0;
400 dnisson 1.8 for (int i = 1; i <= hist_sub->GetXaxis()->GetNbins(); i++) {
401     for (int j = 1; j <= hist_sub->GetYaxis()->GetNbins(); j++) {
402 dnisson 1.9 double nu = energy.bins[i-1][j-1] - sub_energy[i-1][j-1];
403     if (hist_sub->GetBinContent(i, j)
404 dnisson 1.10 - 3.0*pow(chisquare_error(nu), 2.0) > max_sub) {
405 dnisson 1.8 max_sub = hist_sub->GetBinContent(i, j);
406 dnisson 1.9 max_x = static_cast<double>(i-1)*XbinSize
407     + hist_sub->GetXaxis()->GetXmin();
408     max_y = static_cast<double>(j-1)*YbinSize
409     + hist_sub->GetYaxis()->GetXmin();
410 dnisson 1.8 }
411     }
412 dnisson 1.9 }
413 dnisson 1.1
414     for (int i = npar1; i < npar; i++) {
415 dnisson 1.10 if (get_indiv_max_E() == i - npar1) {
416 dnisson 1.9 pval[i] = max_sub / XbinSize / YbinSize;
417 dnisson 1.4 perr[i] = mdef->chisquare_error(pval[i])*0.5;
418 dnisson 1.1 plo[i] = 0.0;
419     phi[i] = 1.0e6;
420     }
421 dnisson 1.10 else if (get_indiv_max_x() == i - npar1) {
422 dnisson 1.1 pval[i] = max_x;
423     perr[i] = (energy.Xhi - energy.Xlo) / static_cast<double>(hist->GetXaxis()->GetNbins());
424     plo[i] = 0.0;
425     phi[i] = 0.0;
426     }
427 dnisson 1.10 else if (get_indiv_max_y() == i - npar1) {
428 dnisson 1.1 pval[i] = max_y;
429     perr[i] = (energy.Yhi - energy.Ylo) / static_cast<double>(hist->GetYaxis()->GetNbins());
430     plo[i] = 0.0;
431     phi[i] = 0.0;
432     }
433     else {
434     string dummy;
435 dnisson 1.10 get_indiv_par(i-npar1, dummy, pval[i], perr[i], plo[i], phi[i]);
436 dnisson 1.1 }
437     }
438     }
439 dnisson 1.10 int n_pars_to_set = init_spec_pars ? ngauss*get_formula()->GetNpar()
440 dnisson 1.1 : npar - npar1;
441     int init_par_to_set = init_spec_pars ? 0 : npar1;
442     for (int i = 0; i < n_pars_to_set; i++) {
443     double dummy;
444     string s;
445     int ipar = i + init_par_to_set;
446 dnisson 1.10 get_indiv_par(i % get_formula()->GetNpar(), s,
447 dnisson 1.1 dummy, dummy, dummy, dummy);
448     ostringstream par_oss;
449 dnisson 1.10 par_oss << s << ipar / (get_formula()->GetNpar()) << flush;
450 dnisson 1.1 try {
451     pars.at(ipar) = TString(par_oss.str());
452     }
453     catch (exception ex) {
454     cerr << "exception 2 in par_init" << endl;
455     }
456     int error_flag;
457     try {
458     gMinuit->mnparm(ipar, pars.at(ipar),
459     pval[ipar], perr[ipar], plo[ipar], phi[ipar], error_flag);
460     } catch (exception ex) {
461     cerr << "exception 3 in par_init" << endl;
462     }
463     }
464     }
465    
466 dnisson 1.11 results fit_histo(TH2 *hist, vector<trouble> &t_vec,
467 dnisson 1.10 void (*cc_minuit)(TMinuit *, TH2 *, int),
468     int start_ngauss,
469     int rebinX, int rebinY,
470     double P_cutoff_val) {
471 dnisson 1.1 TMinuit *gMinuit = new TMinuit();
472     int npar_indiv = mdef->get_formula()->GetNpar();
473     int istat, erflg;
474     results r;
475    
476     gMinuit->SetFCN(fcn);
477     gMinuit->mninit(5, 6, 7);
478    
479     // set error def and machine accuracy
480     double cs_err_def = 1.0;
481     double fcn_eps = 0.2;
482     gMinuit->mnexcm("SET ERR", &cs_err_def, 1, erflg);
483     gMinuit->mnexcm("SET EPS", &fcn_eps, 1, erflg);
484    
485     // rebin histogram
486     TH2 *hist_rebinned;
487     if (rebinX > 1 && rebinY > 1)
488     hist_rebinned = hist->Rebin2D(rebinX, rebinY);
489     else
490     hist_rebinned = hist;
491    
492     r.hist_rebinned = hist_rebinned;
493     // load histogram into energies vector
494     energy.bins.resize(hist_rebinned->GetXaxis()->GetNbins());
495     for (int i = 0; i < energy.bins.size(); i++) {
496     energy.bins[i].resize(hist_rebinned->GetYaxis()->GetNbins());
497     for (int j = 0; j < energy.bins[i].size(); j++) {
498     energy.bins[i][j] = hist_rebinned->GetBinContent(i+1, j+1);
499     }
500     }
501    
502     energy.Xlo = hist->GetXaxis()->GetXmin();
503     energy.Xhi = hist->GetXaxis()->GetXmax();
504     energy.Ylo = hist->GetYaxis()->GetXmin();
505     energy.Yhi = hist->GetYaxis()->GetXmax();
506    
507     if (start_ngauss < 1) {
508     start_ngauss = 1;
509     }
510    
511 dnisson 1.11 for (mdef->set_ngauss(start_ngauss);
512     mdef->get_ngauss() <=
513     (MAX_GAUSS > start_ngauss ? MAX_GAUSS : start_ngauss);
514     mdef->set_ngauss(mdef->get_ngauss()+1)) {
515    
516     int ngauss = mdef->get_ngauss();
517 dnisson 1.1 t_vec.resize(t_vec.size() + 1);
518 dnisson 1.11 int npar = npar_indiv*mdef->get_ngauss();
519 dnisson 1.1 double pval[256], perr[256], plo[256], phi[256];
520     if (npar > 256) {
521     cerr << "Parameter overload" << endl;
522     return r;
523     }
524     vector<TString> pars(npar);
525    
526     // initialize parameters
527 dnisson 1.11 mdef->par_init(hist, gMinuit, pars, pval, perr, plo, phi, npar, r);
528 dnisson 1.5
529 dnisson 1.1 // minimize
530     double chisquare, edm, errdef;
531     int nvpar, nparx;
532     trouble t;
533     t.occ = T_NULL;
534     t.istat = 3;
535 dnisson 1.4 // fix sigmas of the Gaussians
536 dnisson 1.1 for (int i = 0; i < ngauss; i++) {
537     double parno[2];
538 dnisson 1.4 parno[0] = i*npar_indiv + 4;
539     gMinuit->mnexcm("FIX", parno, 1, erflg);
540 dnisson 1.1 }
541     gMinuit->mnexcm("SIMPLEX", 0, 0, erflg);
542     gMinuit->mnexcm("MIGRAD", 0, 0, erflg);
543     gMinuit->mnstat(chisquare, edm, errdef, nvpar, nparx, istat);
544     if (istat == 3) {
545     /* we're not concerned about MINOS errors right now
546     gMinuit->mnexcm("MINOS", 0, 0, erflg);
547     gMinuit->mnstat(chisquare, edm, errdef, nvpar, nparx, istat);
548     if (istat != 3) {
549     t.occ = T_MINOS;
550     t.istat = istat;
551     }
552     */
553     }
554     else {
555     t.occ = T_MIGRAD;
556     t.istat = istat;
557     }
558    
559     try {
560     if (t_vec.size() < ngauss) {
561     t_vec.resize(ngauss);
562     }
563     t_vec.at(ngauss-1) = t;
564     }
565     catch (exception ex) {
566     cerr << "exception in fit_histo" << endl;
567     }
568    
569     // put parameters in map
570     for (int i = 0; i < npar && i < 256; i++) {
571     int iuint; // internal parameter number
572     gMinuit->mnpout(i, pars[i], pval[i], perr[i], plo[i], phi[i], iuint);
573     }
574     vector<double> pval_copy(npar);
575     vector<double> perr_copy(npar);
576     vector<double> plo_copy(npar);
577     vector<double> phi_copy(npar);
578     for (int i = 0; i < npar && i < 256; i++) {
579     pval_copy[i] = pval[i];
580     perr_copy[i] = perr[i];
581     plo_copy[i] = plo[i];
582     phi_copy[i] = phi[i];
583     }
584     r.pars.push_back(pars);
585     r.pval.push_back(pval_copy);
586     r.perr.push_back(perr_copy);
587     r.plo.push_back(plo_copy);
588     r.phi.push_back(phi_copy);
589    
590     // execute user minuit code
591     if (cc_minuit != 0)
592     (*cc_minuit)(gMinuit, hist_rebinned, ngauss);
593    
594     int ndof = 0;
595     if (!ignorezero)
596     ndof = energy.bins.size() * energy.bins[0].size();
597     else {
598     for (int i = 0; i < energy.bins.size(); i++) {
599     for (int j = 0; j < energy.bins[i].size(); j++) {
600     if (energy.bins[i][j] > 1.0e-30)
601     ndof++;
602     }
603     }
604     }
605     ndof -= npar;
606    
607 dnisson 1.2 r.chisquare.push_back(chisquare);
608 dnisson 1.12 double P = TMath::Prob(chisquare, ndof);
609 dnisson 1.1 if (P > P_cutoff_val) {
610     break;
611     }
612     }
613    
614     delete gMinuit;
615     return r;
616     }
617     }