source: trunk/MagicSoft/Mars/mtools/MHSimulatedAnnealing.cc@ 9492

Last change on this file since 9492 was 7749, checked in by tbretz, 19 years ago
*** empty log message ***
File size: 6.1 KB
Line 
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): Markus Gaug <mailto:markus@ifae.es>
19!
20! Copyright: MAGIC Software Development, 2000-2002
21!
22!
23\* ======================================================================== */
24
25///////////////////////////////////////////////////////////////////////
26//
27// MHSimulatedAnnealing
28//
29// This class contains different histograms of the Simulated Annealing
30// Snapshort during optimization and the final results
31//
32///////////////////////////////////////////////////////////////////////
33#include "MHSimulatedAnnealing.h"
34
35#include <TObjArray.h>
36
37#include <TStyle.h>
38#include <TCanvas.h>
39
40#include "MBinning.h"
41
42ClassImp(MHSimulatedAnnealing);
43
44// --------------------------------------------------------------------------
45//
46// Setup histograms
47//
48MHSimulatedAnnealing::MHSimulatedAnnealing(UShort_t moves, UShort_t ndim,
49 const char *name,
50 const char *title)
51 : fDim(ndim), fMoves(moves), fTimeEvolution(NULL), fBestEver(), fBestFuncEval()
52{
53
54 //
55 // set the name and title of this object
56 //
57 fName = name ? name : "MHSimulatedAnnealing";
58 fTitle = title ? title : "Output Histograms of a Simulated Annealing Run";
59
60 fBestEver.SetName("Hist_BestEver");
61 fBestEver.SetTitle("Best Param. Combinations");
62 fBestEver.SetXTitle("Run Duration");
63 fBestEver.SetYTitle("Parameter Nr.");
64 fBestEver.SetZTitle("Parameter Value");
65 fBestEver.SetDirectory(NULL);
66
67 fBestFuncEval.SetName("Hist_BestFuncEval");
68 fBestFuncEval.SetTitle("Best Function Evaluation");
69 fBestFuncEval.SetXTitle("Run Duration");
70 fBestFuncEval.SetYTitle("Function Value");
71 fBestFuncEval.SetDirectory(NULL);
72
73 MBinning binsx, binsy;
74 binsx.SetEdges(fMoves+1, 0, 1);
75 binsy.SetEdges(fDim, 0.5, fDim+0.5);
76 MH::SetBinning(&fBestEver, &binsx, &binsy);
77
78 // For better visibility, omit the first entry in fBestFuncEval
79 // It has nothing significant, anyway
80 binsx.SetEdges(fMoves,1./(fMoves+1), 1);
81 binsx.Apply(fBestFuncEval);
82
83}
84
85void MHSimulatedAnnealing::InitFullSimplex()
86{
87 if (fTimeEvolution)
88 delete fTimeEvolution;
89
90 fTimeEvolution = new TObjArray;
91 fTimeEvolution->SetOwner();
92
93 for (Int_t i=0;i<fDim;i++)
94 {
95 TH2F *hist = new TH2F(Form("Hist_%d", i), Form("Parameter %d", i),
96 fMoves+1, 0., 1.,fDim+1,0.5,fDim+1.5);
97 hist->SetXTitle("Run Duration");
98 hist->SetYTitle("Point Nr. Simplex");
99 hist->SetZTitle(Form("Value of Parameter %d",i));
100 fTimeEvolution->Add(hist);
101 }
102}
103
104Bool_t MHSimulatedAnnealing::StoreFullSimplex(const TMatrix &p, const UShort_t move)
105{
106
107 if (!fTimeEvolution)
108 return kFALSE;
109
110 Int_t idx=0;
111 const Axis_t bin = (move-0.5)/(fMoves+1);
112
113 TIter Next(fTimeEvolution);
114 TH2F *hist=NULL;
115 while ((hist=(TH2F*)Next()))
116 {
117 for (Int_t i=0;i<fDim+1;i++)
118 hist->Fill(bin,i,p(i,idx));
119 idx++;
120 }
121 return kTRUE;
122}
123
124Bool_t MHSimulatedAnnealing::StoreBestValueEver(const TVector &y, const Float_t yb, const UShort_t move)
125{
126 if (y.GetNrows() != fDim)
127 return kFALSE;
128
129 const Axis_t bin = (move-0.5)/(fMoves+1);
130
131 for (Int_t i=0;i<fDim;i++)
132 fBestEver.Fill(bin,0.5+i,((TVector)y)(i));
133
134 fBestFuncEval.Fill(bin,yb);
135
136 return kTRUE;
137}
138
139Bool_t MHSimulatedAnnealing::ChangeTitle(const UShort_t index, const char* title)
140{
141 if (!fTimeEvolution)
142 return kFALSE;
143
144 TH2F *hist = NULL;
145 if (!(hist = (TH2F*)fTimeEvolution->At(index)))
146 return kFALSE;
147
148 hist->SetNameTitle(Form("Hist_%s",title),title);
149 hist->SetYTitle(Form("Value of Parameter %s",title));
150
151 return kTRUE;
152}
153
154void MHSimulatedAnnealing::ChangeFuncTitle(const char* title)
155{
156 fBestFuncEval.SetTitle(title);
157}
158
159TObject *MHSimulatedAnnealing::DrawClone(Option_t *opt) const
160{
161 UShort_t i=2;
162
163 TCanvas *c = MakeDefCanvas(this, 720, 810);
164 if (fTimeEvolution)
165 c->Divide(2,(int)(fDim/2.)+1);
166 else
167 gPad->Divide(1,2);
168
169 gROOT->SetSelectedPad(NULL);
170
171 c->cd(1);
172 gStyle->SetOptStat(0);
173 ((TH1&)fBestFuncEval).DrawCopy();
174
175 c->cd(2);
176 gStyle->SetOptStat(10);
177 ((TH2&)fBestEver).DrawCopy(opt);
178
179 if (fTimeEvolution)
180 {
181 TH2F *hist = NULL;
182 TIter Next(fTimeEvolution);
183 while ((hist=(TH2F*)Next()))
184 {
185 c->cd(++i);
186 hist->DrawCopy(opt);
187 }
188 }
189 c->Modified();
190 c->Update();
191
192 return c;
193}
194
195
196
197// --------------------------------------------------------------------------
198//
199// Draw all histograms.
200//
201void MHSimulatedAnnealing::Draw(Option_t *opt)
202{
203 UShort_t i=2;
204
205 if (!gPad)
206 MakeDefCanvas(this,780,940);
207
208 if (fTimeEvolution)
209 gPad->Divide(2,(int)(fDim/2.)+1);
210 else
211 gPad->Divide(1,2);
212
213 gPad->cd(1);
214 gStyle->SetOptStat(0);
215 fBestFuncEval.Draw();
216 gPad->Modified();
217 gPad->Update();
218
219 gPad->cd(2);
220 gStyle->SetOptStat(10);
221 fBestEver.Draw(opt);
222 gPad->Modified();
223 gPad->Update();
224
225 if (!fTimeEvolution)
226 return;
227
228 TH2F *hist = NULL;
229 TIter Next(fTimeEvolution);
230 while ((hist=(TH2F*)Next()))
231 {
232 gPad->cd(++i);
233 hist->Draw(opt);
234 }
235 gPad->Modified();
236 gPad->Update();
237}
238
239// --------------------------------------------------------------------------
240//
241// Delete the histograms.
242//
243MHSimulatedAnnealing::~MHSimulatedAnnealing()
244{
245 if (fTimeEvolution)
246 delete fTimeEvolution;
247}
248
Note: See TracBrowser for help on using the repository browser.