1 | /* ======================================================================== *\
|
---|
2 | !
|
---|
3 | ! *
|
---|
4 | ! * This file is part of MARS, the MAGIC Analysis and Reconstruction
|
---|
5 | ! * Software. It is distributed to you in the hope that it can be a useful
|
---|
6 | ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
|
---|
7 | ! * It is distributed WITHOUT ANY WARRANTY.
|
---|
8 | ! *
|
---|
9 | ! * Permission to use, copy, modify and distribute this software and its
|
---|
10 | ! * documentation for any purpose is hereby granted without fee,
|
---|
11 | ! * provided that the above copyright notice appear in all copies and
|
---|
12 | ! * that both that copyright notice and this permission notice appear
|
---|
13 | ! * in supporting documentation. It is provided "as is" without express
|
---|
14 | ! * or implied warranty.
|
---|
15 | ! *
|
---|
16 | !
|
---|
17 | !
|
---|
18 | ! Author(s): Thomas Bretz, 3/2004 <mailto:tbretz@astro.uni-wuerzburg.de>
|
---|
19 | !
|
---|
20 | ! Copyright: MAGIC Software Development, 2000-2004
|
---|
21 | !
|
---|
22 | !
|
---|
23 | \* ======================================================================== */
|
---|
24 |
|
---|
25 | /////////////////////////////////////////////////////////////////////////////
|
---|
26 | //
|
---|
27 | // MTFitLoop
|
---|
28 | //
|
---|
29 | // This class optimized parameters which are calculated in an eventloop.
|
---|
30 | // For this it is unimportant whether the loop reads from a file a
|
---|
31 | // matrix or from somewhere else...
|
---|
32 | //
|
---|
33 | // The parameters which are optimized must be stored in a MParContainer
|
---|
34 | // which overwrites SetVariables(). Eg. a MDataPhrase or MF overwrites
|
---|
35 | // SetVariables(). In a MF all arguments given as [0], [1] are
|
---|
36 | // set by SetVariables (in MDataValue).
|
---|
37 | // eg: In you loop you have a cut like this:
|
---|
38 | // MF filter("MHillas.fWidth<[0]");
|
---|
39 | // filter.SetName("MyParameters");
|
---|
40 | //
|
---|
41 | // Now for each time the eventloop is executed
|
---|
42 | // (MEvtLoop::Eventloop(fNumEvents)) the parameters from TMinuit are
|
---|
43 | // passed to MF::SetVariables and changed.
|
---|
44 | //
|
---|
45 | //
|
---|
46 | /////////////////////////////////////////////////////////////////////////////
|
---|
47 | #include "MTFitLoop.h"
|
---|
48 |
|
---|
49 | #include <TArrayD.h>
|
---|
50 | #include <TMinuit.h>
|
---|
51 | #include <TStopwatch.h>
|
---|
52 |
|
---|
53 | #include "MParList.h"
|
---|
54 | #include "MTaskList.h"
|
---|
55 | #include "MEvtLoop.h"
|
---|
56 |
|
---|
57 | #include "MParameters.h"
|
---|
58 |
|
---|
59 | #include "MRead.h"
|
---|
60 |
|
---|
61 | #include "MLog.h"
|
---|
62 | #include "MLogManip.h"
|
---|
63 |
|
---|
64 | ClassImp(MTFitLoop);
|
---|
65 |
|
---|
66 | using namespace std;
|
---|
67 |
|
---|
68 | //------------------------------------------------------------------------
|
---|
69 | //
|
---|
70 | // fcn calculates the function to be minimized (using TMinuit::Migrad)
|
---|
71 | //
|
---|
72 | void MTFitLoop::fcn(Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag)
|
---|
73 | {
|
---|
74 | MTFitLoop *optim = (MTFitLoop*)gMinuit->GetObjectFit();
|
---|
75 |
|
---|
76 | TMinuit *minuit = gMinuit;
|
---|
77 | f = optim->Fcn(npar, gin, par, iflag);
|
---|
78 | gMinuit = minuit;
|
---|
79 |
|
---|
80 | }
|
---|
81 | Double_t MTFitLoop::Fcn(Int_t &npar, Double_t *gin, Double_t *par, Int_t iflag)
|
---|
82 | {
|
---|
83 | MParList *plist = fEvtLoop->GetParList();
|
---|
84 |
|
---|
85 | MParameterD *eval = (MParameterD*) plist->FindCreateObj("MParameterD", fFitParameter);
|
---|
86 | MParContainer *pars = (MParContainer*)plist->FindObject(fParametersName);
|
---|
87 |
|
---|
88 | MRead *read = (MRead*)plist->FindObject("MTaskList")->FindObject("MRead");
|
---|
89 | if (read)
|
---|
90 | read->Rewind();
|
---|
91 |
|
---|
92 | if (fDebug>=1)
|
---|
93 | {
|
---|
94 | Double_t fmin, fedm, errdef;
|
---|
95 | Int_t n1, n2, istat;
|
---|
96 | gMinuit->mnstat(fmin, fedm, errdef, n1, n2, istat);
|
---|
97 | *fLog << inf << underline << "Minimization Status so far:" << endl;
|
---|
98 | *fLog << " Calls: " << gMinuit->fNfcn << endl;
|
---|
99 | *fLog << " Func min: " << fmin << endl;
|
---|
100 | *fLog << " Found edm: " << fedm << endl;
|
---|
101 | *fLog << " ErrDef: " << errdef << endl;
|
---|
102 | *fLog << " Status: ";
|
---|
103 | switch (istat)
|
---|
104 | {
|
---|
105 | case 0: *fLog << "n/a" << endl; break;
|
---|
106 | case 1: *fLog << "approximation only, not accurate" << endl; break;
|
---|
107 | case 2: *fLog << "full matrix, but forced positive-definite" << endl; break;
|
---|
108 | case 3: *fLog << "full accurate covariance matrix" << endl; break;
|
---|
109 | default: *fLog << "undefined" << endl; break;
|
---|
110 | }
|
---|
111 | }
|
---|
112 |
|
---|
113 | if (fDebug>=0)
|
---|
114 | {
|
---|
115 | *fLog << inf << "Set(" << gMinuit->fMaxpar << "): ";
|
---|
116 | for (Int_t i=0; i<gMinuit->fMaxpar; i++)
|
---|
117 | *fLog << par[i] << " ";
|
---|
118 | *fLog << endl;
|
---|
119 | }
|
---|
120 |
|
---|
121 | pars->SetVariables(TArrayD(gMinuit->fMaxpar, par));
|
---|
122 |
|
---|
123 | if (fDebug<3)
|
---|
124 | gLog.SetNullOutput(kTRUE);
|
---|
125 | fEvtLoop->Eventloop(fNumEvents);
|
---|
126 | if (fDebug<3)
|
---|
127 | gLog.SetNullOutput(kFALSE);
|
---|
128 |
|
---|
129 | const Double_t f = eval->GetVal();
|
---|
130 |
|
---|
131 | if (fDebug>=0)
|
---|
132 | *fLog << inf << "F=" << f << endl;
|
---|
133 |
|
---|
134 | if (fDebug>=1)
|
---|
135 | fEvtLoop->GetTaskList()->PrintStatistics();
|
---|
136 |
|
---|
137 | return f;
|
---|
138 | }
|
---|
139 |
|
---|
140 | MTFitLoop::MTFitLoop(Int_t num) : fNum(num), fMaxIterations(1000)
|
---|
141 | {
|
---|
142 | fDebug = -1;
|
---|
143 | fNumEvents = -1;
|
---|
144 | }
|
---|
145 |
|
---|
146 | void MTFitLoop::Optimize(MEvtLoop &loop, TArrayD &pars)
|
---|
147 | {
|
---|
148 | *fLog << inf << "Event loop was setup" << endl;
|
---|
149 | MParList *parlist = loop.GetParList();
|
---|
150 | if (!parlist)
|
---|
151 | return;
|
---|
152 |
|
---|
153 | // MParContainer *pars = (MParContainer*)parlist->FindObject(fParametersName);
|
---|
154 | // if (!pars)
|
---|
155 | // return;
|
---|
156 |
|
---|
157 | fEvtLoop = &loop;
|
---|
158 |
|
---|
159 | // MParContainer ¶meters = *pars;
|
---|
160 |
|
---|
161 | TMinuit *minsave = gMinuit;
|
---|
162 |
|
---|
163 | gMinuit = new TMinuit(pars.GetSize());
|
---|
164 | gMinuit->SetPrintLevel(-1);
|
---|
165 | gMinuit->SetMaxIterations(fMaxIterations);
|
---|
166 |
|
---|
167 | gMinuit->SetFCN(fcn);
|
---|
168 | gMinuit->SetObjectFit(this);
|
---|
169 |
|
---|
170 | // For chisq fits call this // seems to be something like %)
|
---|
171 | //
|
---|
172 | // The default tolerance is 0.1, and the minimization will stop");
|
---|
173 | // when the estimated vertical distance to the minimum (EDM) is");
|
---|
174 | // less than 0.001*[tolerance]*UP (see [SET ERRordef]).");
|
---|
175 | //
|
---|
176 | if (gMinuit->SetErrorDef(1))
|
---|
177 | {
|
---|
178 | *fLog << err << dbginf << "SetErrorDef failed." << endl;
|
---|
179 | return;
|
---|
180 | }
|
---|
181 |
|
---|
182 | //
|
---|
183 | // Set starting values and step sizes for parameters
|
---|
184 | //
|
---|
185 | for (Int_t i=0; i<pars.GetSize(); i++)
|
---|
186 | {
|
---|
187 | TString name = "par[";
|
---|
188 | name += i;
|
---|
189 | name += "]";
|
---|
190 | Double_t vinit = pars[i];
|
---|
191 | Double_t step = fabs(pars[i]/3);
|
---|
192 |
|
---|
193 | Double_t limlo = 0; // limlo=limup=0: no limits
|
---|
194 | Double_t limup = 2*vinit;
|
---|
195 |
|
---|
196 | Bool_t rc = gMinuit->DefineParameter(i, name, vinit, step, 0, limup);
|
---|
197 | if (!rc)
|
---|
198 | continue;
|
---|
199 |
|
---|
200 | *fLog << err << dbginf << "Error in defining parameter #" << i << endl;
|
---|
201 | return;
|
---|
202 | }
|
---|
203 |
|
---|
204 | for (int i=0; i<pars.GetSize() && i<fFixedParams.GetSize(); i++)
|
---|
205 | if (fFixedParams[i]!=0)
|
---|
206 | gMinuit->FixParameter(i);
|
---|
207 |
|
---|
208 | // Now ready for minimization step:
|
---|
209 |
|
---|
210 | TStopwatch clock;
|
---|
211 | clock.Start();
|
---|
212 | const Bool_t rc = gMinuit->Migrad();
|
---|
213 | clock.Stop();
|
---|
214 | clock.Print();
|
---|
215 |
|
---|
216 | if (rc)
|
---|
217 | {
|
---|
218 | *fLog << err << dbginf << "Migrad failed." << endl;
|
---|
219 | return;
|
---|
220 | }
|
---|
221 |
|
---|
222 | *fLog << inf << "Resulting Chisq: " << gMinuit->fAmin << endl;
|
---|
223 |
|
---|
224 | //
|
---|
225 | // Update values of fA, fB:
|
---|
226 | //
|
---|
227 | for (Int_t i=0; i<pars.GetSize(); i++)
|
---|
228 | {
|
---|
229 | Double_t x1, x2;
|
---|
230 | gMinuit->GetParameter(i,x1,x2);
|
---|
231 | pars[i] = x1;
|
---|
232 | cout << i << ": " << pars[i] << endl;
|
---|
233 | }
|
---|
234 |
|
---|
235 | //list.SetVariables(pars);
|
---|
236 |
|
---|
237 | gMinuit = minsave;
|
---|
238 | }
|
---|