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

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