source: trunk/MagicSoft/Mars/mtemp/mmpi/MUnfoldSpectrum.cc@ 6437

Last change on this file since 6437 was 6384, checked in by rwagner, 20 years ago
*** empty log message ***
File size: 9.3 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 expressed
14! * or implied warranty.
15! *
16!
17! Author(s) : R. Wagner, 02/2004 <mailto:rwagner@mppmu.mpg.de>
18!
19! Copyright: MAGIC Software Development, 2000-2005
20!
21!
22\* ======================================================================== */
23
24/////////////////////////////////////////////////////////////////////////////
25//
26// MUnfoldSpectrum
27//
28/////////////////////////////////////////////////////////////////////////////
29#include "MUnfoldSpectrum.h"
30
31#include "TH1D.h"
32#include "TH2D.h"
33#include "TH3D.h"
34
35#include "MLog.h"
36#include "MLogManip.h"
37#include "MUnfold.h"
38
39ClassImp(MUnfoldSpectrum);
40using namespace std;
41
42MUnfoldSpectrum::MUnfoldSpectrum()
43 : fUnfoldingMethod(2), fPrior(0)
44{
45
46}
47
48void MUnfoldSpectrum::Calc()
49{
50 // Unfold # Excess events vs. Energy and Theta
51
52 // TH2D* tobeunfolded = hex->GetHist();
53 // TH2D* unfolded = new TH2D(*tobeunfolded);
54 // TH3D* migration = migm->GetHist();
55
56 const Int_t energyBins = fToBeUnfolded->GetXaxis()->GetNbins();
57 const Int_t thetaBins = fToBeUnfolded->GetYaxis()->GetNbins();
58 const TAxis* axisEnergy = fToBeUnfolded->GetXaxis();
59 const TAxis* axisTheta = fToBeUnfolded->GetYaxis();
60 cout << "Distribution to be unfolded has " << energyBins << ", " << thetaBins << " bins" <<endl;
61
62 TAxis &taxis = *fToBeUnfolded->GetYaxis();
63 Int_t numybins = taxis.GetNbins();
64
65 cout << "Processing a total number of " << numybins << " bins. " <<endl;
66 cout << "-------------------" << endl << endl;
67
68 for (Int_t m=1; m<=numybins; m++) {
69 TString bintitle = "Bin ";
70 bintitle += m;
71 bintitle += ": ";
72
73 cout << "Processing " << bintitle << endl;
74
75 // -----------------------------------------
76 // ha : distribution to be unfolded
77
78 TH1D &ha = *fToBeUnfolded->ProjectionX("", m, m, "e");
79 TString title = bintitle;
80 title += "E-est distr. to be unfolded";
81 ha.SetNameTitle("ha", title);
82 TAxis &aaxis = *ha.GetXaxis();
83 Int_t na = aaxis.GetNbins();
84 Double_t alow = aaxis.GetBinLowEdge(1);
85 Double_t aup = aaxis.GetBinLowEdge(na+1);
86
87 cout << ha.GetName() << ": " << ha.GetTitle() << endl;
88 cout << "-----------------------------------------------------" << endl;
89 for (Int_t i=1; i<=ha.GetNbinsX(); i++)
90 cout << ha.GetBinContent(i) << " \t";
91 cout << endl << endl;
92
93 // -----------------------------------------
94 // covariance matrix of the distribution ha
95
96 title = bintitle;
97 title += "Error matrix of distribution ha";
98 TH2D hacov("hacov", title, na, alow, aup, na, alow, aup);
99
100 Double_t errmin = 3.0;
101 for (Int_t i=1; i<=na; i++) {
102 for (Int_t j=1; j<=na; j++)
103 hacov.SetBinContent(i, j, 0.0);
104 const Double_t content = ha.GetBinContent(i);
105 const Double_t error2 = (ha.GetBinError(i))*(ha.GetBinError(i));
106 if (content <= errmin && error2 < errmin)
107 hacov.SetBinContent(i, i, errmin);
108 else
109 hacov.SetBinContent(i, i, error2);
110 }
111
112 //PrintTH2Content(hacov);
113
114 // -----------------------------------------
115 // migration matrix :
116 // x corresponds to measured quantity
117 // y corresponds to true quantity
118
119 // The projection is made for the selected bins only.
120 // To select a bin range along an axis, use TAxis::SetRange, eg
121 // h3.GetYaxis()->SetRange(23,56);
122
123 // taxis->SetRange(m,m);
124 TH2D &hmig = *(TH2D*)fMigrationMatrix->Project3D("yxe");
125 title = bintitle;
126 title += "Migration Matrix";
127 hmig.SetNameTitle("Migrat", title);
128
129 TAxis &aaxismig = *hmig.GetXaxis();
130 Int_t namig = aaxismig.GetNbins();
131
132 if (na != namig) {
133 cout << "doUnfolding : binnings are incompatible; na, namig = "
134 << na << ", " << namig << endl;
135 return;
136 }
137
138 TAxis &baxismig = *hmig.GetYaxis();
139 Int_t nbmig = baxismig.GetNbins();
140 Double_t blow = baxismig.GetBinLowEdge(1);
141 Double_t bup = baxismig.GetBinLowEdge(nbmig+1);
142
143
144 // -----------------------------------------
145 // dummy ideal distribution
146
147 Int_t nb = nbmig;
148
149 title = bintitle;
150 title += "Dummy Ideal distribution";
151 TH1D hb0("dummyhb0", title, nb, blow, bup);
152 //MH::SetBinning(&hb0, &baxismig);
153 hb0.Sumw2();
154
155 for (Int_t k=1; k<=nb; k++) {
156 hb0.SetBinContent(k, 1.0/nb);
157 hb0.SetBinError (k, 0.1/nb);
158 }
159
160 //PrintTH1Content(hb0);
161
162 // -----------------------------------------
163 // unfolded distribution
164
165 title = bintitle;
166 title += "Unfolded distribution";
167 TH1D hb("hb", title, nb, blow, bup);
168
169 // ha is the distribution to be unfolded
170 // hacov is the covariance matrix of the distribution ha
171 // hmig is the migration matrix;
172 // it is used in the unfolding unless it is overwritten
173 // by SmoothMigrationMatrix by the smoothed migration matrix
174 // hmigor is the migration matrix to be smoothed;
175 // the smoothed migration matrix will be used in the unfolding
176 // hpr the prior distribution
177 // it is only used if SetPriorInput(*hpr) is called
178 // hb unfolded distribution
179
180 // create an MUnfold object;
181 MUnfold unfold(ha, hacov, hmig);
182 unfold.bintitle = bintitle;
183
184 // smooth the migration matrix;
185 // the smoothed migration matrix will be used in the unfolding
186 // hmig is the original (unsmoothed) migration matrix
187 unfold.SmoothMigrationMatrix(hmig);
188
189 // define prior distribution (has always to be defined)
190 // the alternatives are
191 // SetPriorConstant(): isotropic distribution
192 // SetPriorPower(gamma): dN/dE = E^{-gamma}
193 // SetPriorInput(*hpr): the distribution *hpr is used
194 // SetPriorRebin(*ha): use rebinned histogram ha
195 Bool_t errorprior=kTRUE;
196 switch (fPrior)
197 {
198 case 1:
199 unfold.SetPriorConstant();
200 break;
201 case 2:
202 errorprior = unfold.SetPriorPower(fPriorPowerGamma);
203 break;
204 case 3:
205 if (!fPriorInputHist)
206 {
207 cout << "Error: No hpr!" << endl;
208 return;
209 }
210 errorprior = unfold.SetPriorInput(*fPriorInputHist);
211 break;
212 case 4:
213 errorprior = unfold.SetPriorRebin(*fPriorRebinHist);
214 break;
215 }
216 if (!errorprior)
217 {
218 cout << "MUnfoldSpectrum::SetPrior... : failed. fPrior = " ;
219 cout << fPrior << endl;
220 return;
221 }
222
223 // calculate the matrix G = M * M(transposed)
224 // M being the migration matrix
225 unfold.CalculateG();
226
227 switch (fUnfoldingMethod)
228 {
229 case 1: // Schmelling:
230 // minimize the function Z by Gauss-Newton iteration;
231 // the parameters to be fitted are gamma(i) = lambda(i)/w;
232 cout << "Unfolding algorithm : Schmelling" << endl;
233 if (!unfold.Schmelling(hb0)) cout << "MUnfoldSpectrum::Schmelling : failed." << endl;
234 break;
235
236 case 2: // Tikhonov2 :
237 // regularization term is sum of (2nd deriv.)**2 ;
238 // minimization by using MINUIT;
239 // the parameters to be fitted are
240 // the bin contents of the unfolded distribution
241 cout << "Unfolding algorithm : Tikhonov" << endl;
242 if (!unfold.Tikhonov2(hb0)) cout << "MUnfoldSpectrum::Tikhonov2 : failed." << endl;
243 break;
244
245 case 3: // Bertero: minimization by iteration
246 cout << "Unfolding algorithm : Bertero" << endl;
247 if (!unfold.Bertero(hb0)) cout << "MUnfoldSpectrum::Bertero : failed." << endl;
248 break;
249 }
250 unfold.PrintResults();
251 unfold.DrawPlots();
252
253 // get unfolded distribution
254 TMatrixD &Vb = unfold.GetVb();
255 TMatrixD &Vbcov = unfold.GetVbcov();
256 UInt_t fNb = unfold.fNb;
257 for (UInt_t a=0; a<fNb; a++)
258 {
259 hb.SetBinContent(a+1, Vb(a,0));
260 hb.SetBinError(a+1, sqrt(Vbcov(a, a)) );
261 }
262
263
264 for (Int_t k=1; k<=nb; k++) {
265 Double_t content = hb.GetBinContent(k);
266 Double_t error = hb.GetBinError(k);
267
268 fUnfolded->SetBinContent(k, m, content);
269 fUnfolded->SetBinError(k, m, error);
270
271 //hex->FillBinContent(k, m, content, error);
272 }
273 }
274
275
276
277
278
279}
280
281
282// --------------------------------------------------------------------------
283//
284// Default destructor.
285//
286MUnfoldSpectrum::~MUnfoldSpectrum()
287{
288}
289
290
291// -----------------------------------------------------------------------
292//
293// Define prior distribution to be a constant
294//
295void MUnfoldSpectrum::SetPriorConstant()
296{
297 fPrior=1;
298}
299
300void MUnfoldSpectrum::SetPriorRebin(TH1D *ha)
301{
302 fPriorRebinHist=ha;
303 fPrior=2;
304}
305
306void MUnfoldSpectrum::SetPriorInput(TH1D *hpr)
307{
308 fPriorInputHist=hpr;
309 fPrior=3;
310}
311
312void MUnfoldSpectrum::SetPriorPower(Double_t gamma)
313{
314 fPriorPowerGamma=gamma;
315 fPrior=4;
316}
317
Note: See TracBrowser for help on using the repository browser.