| 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 | ! Author(s): Javier Rico 02/2005 <mailto:jrico@ifae.es>
|
|---|
| 18 | !
|
|---|
| 19 | ! Copyright: MAGIC Software Development, 2000-2005
|
|---|
| 20 | !
|
|---|
| 21 | !
|
|---|
| 22 | \* ======================================================================== */
|
|---|
| 23 |
|
|---|
| 24 | //////////////////////////////////////////////////////////////////////////////
|
|---|
| 25 | //
|
|---|
| 26 | // MHMcUnfoldCoeff
|
|---|
| 27 | // Container that stores the coefficients to convert from estimated to real
|
|---|
| 28 | // (unfolded) energy.
|
|---|
| 29 | // It is filled by MMcUnfoldCoeffCalc
|
|---|
| 30 | // Events are weighted from any initial spectrum to a tentative spectrum
|
|---|
| 31 | // that must be set through MMcUnfoldCoeffCalc::SetSpectrum method.
|
|---|
| 32 | //
|
|---|
| 33 | //////////////////////////////////////////////////////////////////////////////
|
|---|
| 34 |
|
|---|
| 35 | #include <fstream>
|
|---|
| 36 | #include <math.h>
|
|---|
| 37 |
|
|---|
| 38 | #include "TH1D.h"
|
|---|
| 39 | #include "TH2D.h"
|
|---|
| 40 | #include "TF1.h"
|
|---|
| 41 | #include "TCanvas.h"
|
|---|
| 42 | #include "TArrayD.h"
|
|---|
| 43 |
|
|---|
| 44 | #include "MHMcUnfoldCoeff.h"
|
|---|
| 45 |
|
|---|
| 46 | #include "MBinning.h"
|
|---|
| 47 |
|
|---|
| 48 | #include "MLog.h"
|
|---|
| 49 | #include "MLogManip.h"
|
|---|
| 50 |
|
|---|
| 51 | ClassImp(MHMcUnfoldCoeff);
|
|---|
| 52 |
|
|---|
| 53 | using namespace std;
|
|---|
| 54 |
|
|---|
| 55 |
|
|---|
| 56 | // -------------------------------------------------------------------------
|
|---|
| 57 | //
|
|---|
| 58 | // Constructor: Define and configure the different histograms to be filled
|
|---|
| 59 | //
|
|---|
| 60 | MHMcUnfoldCoeff::MHMcUnfoldCoeff(const char *name, const char *title) :
|
|---|
| 61 | fNsubbins(20), fEmin(1.), fEmax(100000.), fNumMin(10), fFineBinning(NULL)
|
|---|
| 62 | {
|
|---|
| 63 | fName = name ? name : "MHMcUnfoldCoeff";
|
|---|
| 64 | fTitle = title ? title : "Unfolding Coefficients vs. Theta vs. Energy";
|
|---|
| 65 |
|
|---|
| 66 | // define histos
|
|---|
| 67 | fHistAll = new TH1D();
|
|---|
| 68 | fHistWeight = new TH1D();
|
|---|
| 69 | fHistMcE = new TH2D();
|
|---|
| 70 | fHistEstE = new TH2D();
|
|---|
| 71 | fHistCoeff = new TH2D();
|
|---|
| 72 |
|
|---|
| 73 | // set histo names
|
|---|
| 74 | fHistAll->SetName("AllEvents");
|
|---|
| 75 | fHistWeight->SetName("Weights");
|
|---|
| 76 | fHistMcE->SetName("MCEnergy");
|
|---|
| 77 | fHistEstE->SetName("EstimatedEnergy");
|
|---|
| 78 | fHistCoeff->SetName("UnfoldCoeff");
|
|---|
| 79 |
|
|---|
| 80 | // set histo titles
|
|---|
| 81 | fHistAll->SetTitle("MC energy for all generated events");
|
|---|
| 82 | fHistWeight->SetTitle("Weights");
|
|---|
| 83 | fHistMcE->SetTitle("MC energy for survivors");
|
|---|
| 84 | fHistEstE->SetTitle("Estimate energy for survivors");
|
|---|
| 85 | fHistCoeff->SetTitle(fTitle);
|
|---|
| 86 |
|
|---|
| 87 | // histos directory
|
|---|
| 88 | fHistAll->SetDirectory(NULL);
|
|---|
| 89 | fHistWeight->SetDirectory(NULL);
|
|---|
| 90 | fHistMcE->SetDirectory(NULL);
|
|---|
| 91 | fHistEstE->SetDirectory(NULL);
|
|---|
| 92 | fHistCoeff->SetDirectory(NULL);
|
|---|
| 93 |
|
|---|
| 94 | // histos style
|
|---|
| 95 | fHistAll->UseCurrentStyle();
|
|---|
| 96 | fHistWeight->UseCurrentStyle();
|
|---|
| 97 | fHistMcE->UseCurrentStyle();
|
|---|
| 98 | fHistEstE->UseCurrentStyle();
|
|---|
| 99 | fHistCoeff->UseCurrentStyle();
|
|---|
| 100 |
|
|---|
| 101 | // labels
|
|---|
| 102 | fHistAll->SetXTitle("E [GeV]");
|
|---|
| 103 | fHistWeight->SetXTitle("E [GeV]");
|
|---|
| 104 | fHistMcE->SetXTitle("E_{MC} [GeV]");
|
|---|
| 105 | fHistEstE->SetXTitle("E_{est} [GeV]");
|
|---|
| 106 | fHistCoeff->SetXTitle("E_{est} [GeV]");
|
|---|
| 107 |
|
|---|
| 108 | fHistWeight->SetYTitle("weight");
|
|---|
| 109 | fHistMcE->SetYTitle("\\theta [\\circ]");
|
|---|
| 110 | fHistEstE->SetYTitle("\\theta [\\circ]");
|
|---|
| 111 | fHistCoeff->SetYTitle("\\theta [\\circ]");
|
|---|
| 112 |
|
|---|
| 113 | fHistMcE->SetZTitle("weighted entries");
|
|---|
| 114 | fHistEstE->SetZTitle("weighted entries");
|
|---|
| 115 | fHistCoeff->SetZTitle("Coefficient");
|
|---|
| 116 | }
|
|---|
| 117 |
|
|---|
| 118 | // -------------------------------------------------------------------------
|
|---|
| 119 | //
|
|---|
| 120 | // Destructor: Delete all the histograms
|
|---|
| 121 | //
|
|---|
| 122 | MHMcUnfoldCoeff::~MHMcUnfoldCoeff()
|
|---|
| 123 | {
|
|---|
| 124 | delete fHistAll;
|
|---|
| 125 | delete fHistWeight;
|
|---|
| 126 | delete fHistMcE;
|
|---|
| 127 | delete fHistEstE;
|
|---|
| 128 | delete fHistCoeff;
|
|---|
| 129 |
|
|---|
| 130 | if(fFineBinning)
|
|---|
| 131 | delete fFineBinning;
|
|---|
| 132 | }
|
|---|
| 133 |
|
|---|
| 134 | // --------------------------------------------------------------------------
|
|---|
| 135 | //
|
|---|
| 136 | // Set the histograms binning. The coarse binning has to be provided
|
|---|
| 137 | // from outside (in the parameter list) whereas the fine binning for
|
|---|
| 138 | // energy is computed from the former.
|
|---|
| 139 | // The subbinning is set logarithmicaly (linearly) for a logarithmic
|
|---|
| 140 | // (linear) coarse binning.
|
|---|
| 141 | // It is extended down to fEmin and up to fEmax
|
|---|
| 142 | // The number of subbins per coarse bins may be changed using
|
|---|
| 143 | // SetNsubbins() method
|
|---|
| 144 | // The maximum and minimum energies of the fine binning histograms
|
|---|
| 145 | // may be changed with the SetEmin() and SetEmax() methods
|
|---|
| 146 | //
|
|---|
| 147 | void MHMcUnfoldCoeff::SetCoarseBinnings(const MBinning &binsEnergy,
|
|---|
| 148 | const MBinning &binsTheta)
|
|---|
| 149 | {
|
|---|
| 150 | // set histo binning (for coarse bin histos)
|
|---|
| 151 | MH::SetBinning(fHistMcE,&binsEnergy,&binsTheta);
|
|---|
| 152 | MH::SetBinning(fHistEstE,&binsEnergy,&binsTheta);
|
|---|
| 153 | MH::SetBinning(fHistCoeff,&binsEnergy,&binsTheta);
|
|---|
| 154 |
|
|---|
| 155 | // good errors
|
|---|
| 156 | fHistMcE->Sumw2();
|
|---|
| 157 | fHistEstE->Sumw2();
|
|---|
| 158 |
|
|---|
| 159 | // define fine binning from coarse one
|
|---|
| 160 | const Int_t nbins = binsEnergy.GetNumBins();
|
|---|
| 161 | const TArrayD edges = binsEnergy.GetEdgesD();
|
|---|
| 162 | const Bool_t isLinear = binsEnergy.IsLinear();
|
|---|
| 163 |
|
|---|
| 164 | // compute bins needed to extend down to fEmin and up to fEmax
|
|---|
| 165 | Double_t dedown, deup;
|
|---|
| 166 | Int_t binsdown,binsup;
|
|---|
| 167 | MBinning dbin;
|
|---|
| 168 | MBinning ubin;
|
|---|
| 169 | if(isLinear)
|
|---|
| 170 | {
|
|---|
| 171 | dedown = (edges[1]-edges[0])/fNsubbins;
|
|---|
| 172 | deup = (edges[nbins]-edges[nbins-1])/fNsubbins;
|
|---|
| 173 | binsdown = (fEmin<edges[0]) ? Int_t((edges[0]-fEmin)/dedown) : 0;
|
|---|
| 174 | binsup = (fEmax>edges[nbins]) ? Int_t((fEmax-edges[nbins])/deup) : 0;
|
|---|
| 175 | if(binsdown) dbin.SetEdges(binsdown,fEmin,edges[0]);
|
|---|
| 176 | if(binsup) ubin.SetEdges(binsup,edges[nbins],fEmax);
|
|---|
| 177 | }
|
|---|
| 178 | else
|
|---|
| 179 | {
|
|---|
| 180 | dedown = (TMath::Log10(edges[1])-TMath::Log10(edges[0]))/fNsubbins;
|
|---|
| 181 | deup = (TMath::Log10(edges[nbins])-TMath::Log10(edges[nbins-1]))/fNsubbins;
|
|---|
| 182 | binsdown = (fEmin<edges[0]) ? Int_t((TMath::Log10(edges[0])-TMath::Log10(fEmin))/dedown) : 0;
|
|---|
| 183 | binsup = (fEmax>edges[nbins])? Int_t((TMath::Log10(fEmax)-TMath::Log10(edges[nbins]))/deup) : 0;
|
|---|
| 184 | if(binsdown) dbin.SetEdgesLog(binsdown,fEmin,edges[0]);
|
|---|
| 185 | if(binsup) ubin.SetEdgesLog(binsup,edges[nbins],fEmax);
|
|---|
| 186 | }
|
|---|
| 187 |
|
|---|
| 188 | // fill the subins' edges
|
|---|
| 189 | TArrayD fineedges(binsdown+nbins*fNsubbins+binsup+1);
|
|---|
| 190 |
|
|---|
| 191 | for(Int_t i=0;i<binsdown;i++)
|
|---|
| 192 | fineedges[i] = dbin.GetEdgesD().At(i);
|
|---|
| 193 |
|
|---|
| 194 | for(Int_t i=0;i<nbins;i++)
|
|---|
| 195 | {
|
|---|
| 196 | MBinning coarb;
|
|---|
| 197 | if(isLinear)
|
|---|
| 198 | coarb.SetEdges(fNsubbins,edges[i],edges[i+1]);
|
|---|
| 199 | else
|
|---|
| 200 | coarb.SetEdgesLog(fNsubbins,edges[i],edges[i+1]);
|
|---|
| 201 |
|
|---|
| 202 | for(Int_t j=0;j<fNsubbins;j++)
|
|---|
| 203 | fineedges[binsdown+i*fNsubbins+j] = coarb.GetEdgesD().At(j);
|
|---|
| 204 | }
|
|---|
| 205 |
|
|---|
| 206 | for(Int_t i=0;i<binsup;i++)
|
|---|
| 207 | fineedges[binsdown+nbins*fNsubbins+i] = ubin.GetEdgesD().At(i);
|
|---|
| 208 |
|
|---|
| 209 | fineedges[binsdown+nbins*fNsubbins+binsup] = binsup ? ubin.GetEdgesD().At(binsup): edges[nbins];
|
|---|
| 210 |
|
|---|
| 211 | // define fine binning object
|
|---|
| 212 | fFineBinning = new MBinning;
|
|---|
| 213 | fFineBinning->SetEdges(fineedges);
|
|---|
| 214 |
|
|---|
| 215 | // apply fine binning to histograms that need it
|
|---|
| 216 | fFineBinning->Apply(*fHistAll);
|
|---|
| 217 | fFineBinning->Apply(*fHistWeight);
|
|---|
| 218 | }
|
|---|
| 219 |
|
|---|
| 220 | // -------------------------------------------------------------------------
|
|---|
| 221 | //
|
|---|
| 222 | // Compute the weights for a particular tentative spectrum.
|
|---|
| 223 | // For each energy (fine) bin we compute it as the value of the input function
|
|---|
| 224 | // divided by the number of entries in the actual energy histogram.
|
|---|
| 225 | // First and last bin are not filled since they could be biased
|
|---|
| 226 | //
|
|---|
| 227 | void MHMcUnfoldCoeff::ComputeWeights(TF1* spectrum)
|
|---|
| 228 | {
|
|---|
| 229 | Bool_t first=kTRUE;
|
|---|
| 230 | Int_t lastb=0;
|
|---|
| 231 | for(Int_t i=0;i<fHistAll->GetNbinsX();i++)
|
|---|
| 232 | {
|
|---|
| 233 | const Float_t denom = fHistAll->GetBinContent(i+1); // number of events in initial spectrum
|
|---|
| 234 | const Float_t ew = fHistAll->GetBinCenter(i+1); // energy associated to the bin
|
|---|
| 235 | const Float_t numer = spectrum->Eval(ew); // number of events for the required spectrum
|
|---|
| 236 |
|
|---|
| 237 | if(denom>fNumMin)
|
|---|
| 238 | {
|
|---|
| 239 | // do not fill it if it is the first one
|
|---|
| 240 | if(first)
|
|---|
| 241 | {
|
|---|
| 242 | fHistWeight->SetBinContent(i+1,-1);
|
|---|
| 243 | first=kFALSE;
|
|---|
| 244 | }
|
|---|
| 245 | else
|
|---|
| 246 | {
|
|---|
| 247 | fHistWeight->SetBinContent(i+1,numer/denom);
|
|---|
| 248 | lastb=i+1;
|
|---|
| 249 | }
|
|---|
| 250 | }
|
|---|
| 251 | else
|
|---|
| 252 | fHistWeight->SetBinContent(i+1,-1);
|
|---|
| 253 | }
|
|---|
| 254 |
|
|---|
| 255 | //remove last filled bin
|
|---|
| 256 | if(lastb)
|
|---|
| 257 | fHistWeight->SetBinContent(lastb,-1);
|
|---|
| 258 | }
|
|---|
| 259 |
|
|---|
| 260 | // --------------------------------------------------------------
|
|---|
| 261 | //
|
|---|
| 262 | // Compute the coefficients used for the (iterative) unfolding
|
|---|
| 263 | // The coefficients are the ratio between the contents of the
|
|---|
| 264 | // mc energy and estimate energy histograms (filled with weighted
|
|---|
| 265 | // events
|
|---|
| 266 | // Errors are computed as if estimated and MC energy histos where
|
|---|
| 267 | // uncorrelated
|
|---|
| 268 | //
|
|---|
| 269 | void MHMcUnfoldCoeff::ComputeCoefficients()
|
|---|
| 270 | {
|
|---|
| 271 | for(Int_t j=0;j<fHistEstE->GetNbinsY();j++)
|
|---|
| 272 | for(Int_t i=0;i<fHistEstE->GetNbinsX();i++)
|
|---|
| 273 | {
|
|---|
| 274 | const Float_t num = fHistMcE->GetBinContent(i+1,j+1);
|
|---|
| 275 | const Float_t den = fHistEstE->GetBinContent(i+1,j+1);
|
|---|
| 276 |
|
|---|
| 277 | if(den>0)
|
|---|
| 278 | {
|
|---|
| 279 | const Float_t cf = num/den;
|
|---|
| 280 | fHistCoeff->SetBinContent(i+1,j+1,cf);
|
|---|
| 281 |
|
|---|
| 282 | if(num>0)
|
|---|
| 283 | {
|
|---|
| 284 | const Float_t ernum = fHistMcE->GetBinError(i+1,j+1);
|
|---|
| 285 | const Float_t erden = fHistEstE->GetBinError(i+1,j+1);
|
|---|
| 286 | const Float_t ercf = cf*(TMath::Sqrt(ernum/num*ernum/num+erden/den*erden/den));
|
|---|
| 287 | fHistCoeff->SetBinError(i+1,j+1,ercf);
|
|---|
| 288 | }
|
|---|
| 289 | else
|
|---|
| 290 | fHistCoeff->SetBinError(i+1,j+1,0);
|
|---|
| 291 | }
|
|---|
| 292 | else
|
|---|
| 293 | {
|
|---|
| 294 | *fLog << warn << "MHMcUnfoldCoeff::ComputeCoefficients Warning: energy bin " << i << ", thetabin " << j << " has no survivors, setting coefficient to 0" << endl;
|
|---|
| 295 | fHistCoeff->SetBinContent(i+1,j+1,0.);
|
|---|
| 296 | fHistCoeff->SetBinError(i+1,j+1,0.);
|
|---|
| 297 | }
|
|---|
| 298 | }
|
|---|
| 299 | }
|
|---|
| 300 |
|
|---|
| 301 | // --------------------------------------------------------------------------
|
|---|
| 302 | //
|
|---|
| 303 | // Fill data into the histogram which contains all showers
|
|---|
| 304 | //
|
|---|
| 305 | void MHMcUnfoldCoeff::FillAll(Double_t energy)
|
|---|
| 306 | {
|
|---|
| 307 | fHistAll->Fill(energy);
|
|---|
| 308 | }
|
|---|
| 309 |
|
|---|
| 310 | // --------------------------------------------------------------------------
|
|---|
| 311 | //
|
|---|
| 312 | // Fill data into the histograms which contain survivors
|
|---|
| 313 | // each event is introduced with a weight depending on its (MC) energy
|
|---|
| 314 | //
|
|---|
| 315 | void MHMcUnfoldCoeff::FillSel(Double_t mcenergy,Double_t estenergy,Double_t theta)
|
|---|
| 316 | {
|
|---|
| 317 | // bin of the corresponding weight
|
|---|
| 318 | const Int_t wbin = fFineBinning->FindHiEdge(mcenergy);
|
|---|
| 319 |
|
|---|
| 320 | if(wbin>0)
|
|---|
| 321 | {
|
|---|
| 322 | // weight
|
|---|
| 323 | const Double_t weight = fHistWeight->GetBinContent(wbin);
|
|---|
| 324 |
|
|---|
| 325 | if(weight>0)
|
|---|
| 326 | {
|
|---|
| 327 | fHistMcE->Fill(mcenergy,theta,weight);
|
|---|
| 328 | fHistEstE->Fill(estenergy,theta,weight);
|
|---|
| 329 | }
|
|---|
| 330 | // else
|
|---|
| 331 | // *fLog << warn << "MHMcUnfoldCoeff::FillSel Warning: event with energy " << mcenergy << " has no computed weight (lack of MC statistics), skipping" << endl;
|
|---|
| 332 | }
|
|---|
| 333 | // the event has an energy out of the considered range
|
|---|
| 334 | else
|
|---|
| 335 | *fLog << warn << "MHMcUnfoldCoeff::FillSel Warning: event with energy " << mcenergy << " has energy out of bounds, skipping" << endl;
|
|---|
| 336 | }
|
|---|
| 337 |
|
|---|
| 338 | // --------------------------------------------------------------------------
|
|---|
| 339 | //
|
|---|
| 340 | // Draw
|
|---|
| 341 | //
|
|---|
| 342 | void MHMcUnfoldCoeff::Draw(Option_t* option)
|
|---|
| 343 | {
|
|---|
| 344 | TCanvas *c1 = new TCanvas();
|
|---|
| 345 | c1->SetLogx();
|
|---|
| 346 | c1->SetLogy();
|
|---|
| 347 | c1->SetGridx();
|
|---|
| 348 | c1->SetGridy();
|
|---|
| 349 |
|
|---|
| 350 | fHistCoeff->Draw("");
|
|---|
| 351 | }
|
|---|
| 352 |
|
|---|