| 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 Hengstebeck 2/2005 <mailto:hengsteb@physik.hu-berlin.de>
|
|---|
| 19 | ! Author(s): Thomas Bretz 8/2005 <mailto:tbretz@astro.uni-wuerzburg.de>
|
|---|
| 20 | !
|
|---|
| 21 | ! Copyright: MAGIC Software Development, 2000-2005
|
|---|
| 22 | !
|
|---|
| 23 | !
|
|---|
| 24 | \* ======================================================================== */
|
|---|
| 25 |
|
|---|
| 26 | /////////////////////////////////////////////////////////////////////////////
|
|---|
| 27 | //
|
|---|
| 28 | // MRFEnergyEst
|
|---|
| 29 | //
|
|---|
| 30 | //
|
|---|
| 31 | ////////////////////////////////////////////////////////////////////////////
|
|---|
| 32 | #include "MRFEnergyEst.h"
|
|---|
| 33 |
|
|---|
| 34 | #include <TVector.h>
|
|---|
| 35 |
|
|---|
| 36 | #include "MHMatrix.h"
|
|---|
| 37 |
|
|---|
| 38 | #include "MLog.h"
|
|---|
| 39 | #include "MLogManip.h"
|
|---|
| 40 |
|
|---|
| 41 | #include "MData.h"
|
|---|
| 42 | #include "MDataArray.h"
|
|---|
| 43 |
|
|---|
| 44 | #include "MRanForest.h"
|
|---|
| 45 | #include "MParameters.h"
|
|---|
| 46 |
|
|---|
| 47 | #include "MParList.h"
|
|---|
| 48 | #include "MTaskList.h"
|
|---|
| 49 | #include "MEvtLoop.h"
|
|---|
| 50 | #include "MRanForestGrow.h"
|
|---|
| 51 |
|
|---|
| 52 | ClassImp(MRFEnergyEst);
|
|---|
| 53 |
|
|---|
| 54 | using namespace std;
|
|---|
| 55 |
|
|---|
| 56 | const TString MRFEnergyEst::gsDefName = "MRFEnergyEst";
|
|---|
| 57 | const TString MRFEnergyEst::gsDefTitle = "RF for energy estimation";
|
|---|
| 58 |
|
|---|
| 59 | MRFEnergyEst::MRFEnergyEst(const char *name, const char *title)
|
|---|
| 60 | : fDebug(kFALSE), fData(0), fEnergyEst(0),
|
|---|
| 61 | fNumTrees(-1), fNumTry(-1), fNdSize(-1),
|
|---|
| 62 | fTestMatrix(0), fEstimationMode(kMean)
|
|---|
| 63 | {
|
|---|
| 64 | fName = name ? name : gsDefName.Data();
|
|---|
| 65 | fTitle = title ? title : gsDefTitle.Data();
|
|---|
| 66 |
|
|---|
| 67 | // FIXME:
|
|---|
| 68 | fNumTrees = 100; //100
|
|---|
| 69 | fNumTry = 5; //3
|
|---|
| 70 | fNdSize = 0; //1 0 means: in MRanForest estimated best value will be calculated
|
|---|
| 71 | }
|
|---|
| 72 |
|
|---|
| 73 | MRFEnergyEst::~MRFEnergyEst()
|
|---|
| 74 | {
|
|---|
| 75 | fEForests.Delete();
|
|---|
| 76 | }
|
|---|
| 77 |
|
|---|
| 78 | Int_t MRFEnergyEst::Train(const MHMatrix &matrixtrain, const TArrayD &grid, Int_t ver)
|
|---|
| 79 | {
|
|---|
| 80 | gLog.Separator("MRFEnergyEst - Train");
|
|---|
| 81 |
|
|---|
| 82 | if (!matrixtrain.GetColumns())
|
|---|
| 83 | {
|
|---|
| 84 | *fLog << err << "ERROR - MHMatrix does not contain rules... abort." << endl;
|
|---|
| 85 | return kFALSE;
|
|---|
| 86 | }
|
|---|
| 87 |
|
|---|
| 88 | const Int_t ncols = matrixtrain.GetM().GetNcols();
|
|---|
| 89 | const Int_t nrows = matrixtrain.GetM().GetNrows();
|
|---|
| 90 | if (ncols<=0 || nrows <=0)
|
|---|
| 91 | {
|
|---|
| 92 | *fLog << err << "ERROR - No. of columns or no. of rows of matrixtrain equal 0 ... abort." << endl;
|
|---|
| 93 | return kFALSE;
|
|---|
| 94 | }
|
|---|
| 95 |
|
|---|
| 96 | // rules (= combination of image par) to be used for energy estimation
|
|---|
| 97 | TFile fileRF(fFileName, "recreate");
|
|---|
| 98 | if (!fileRF.IsOpen())
|
|---|
| 99 | {
|
|---|
| 100 | *fLog << err << "ERROR - File to store RFs could not be opened... abort." << endl;
|
|---|
| 101 | return kFALSE;
|
|---|
| 102 | }
|
|---|
| 103 |
|
|---|
| 104 | const Int_t nobs = 3; // Number of obsolete columns
|
|---|
| 105 |
|
|---|
| 106 | const MDataArray &dcol = *matrixtrain.GetColumns();
|
|---|
| 107 |
|
|---|
| 108 | MDataArray usedrules;
|
|---|
| 109 | for (Int_t i=0; i<ncols-nobs; i++) // -3 is important!!!
|
|---|
| 110 | usedrules.AddEntry(dcol[i].GetRule());
|
|---|
| 111 |
|
|---|
| 112 | MDataArray rules(usedrules);
|
|---|
| 113 | rules.AddEntry(ver<2?"Classification":dcol[ncols-1].GetRule());
|
|---|
| 114 |
|
|---|
| 115 | // prepare matrix for current energy bin
|
|---|
| 116 | TMatrix mat(matrixtrain.GetM());
|
|---|
| 117 |
|
|---|
| 118 | // last column must be removed (true energy col.)
|
|---|
| 119 | mat.ResizeTo(nrows, ncols-nobs+1);
|
|---|
| 120 |
|
|---|
| 121 | if (fDebug)
|
|---|
| 122 | gLog.SetNullOutput(kTRUE);
|
|---|
| 123 |
|
|---|
| 124 | const Int_t nbins = ver>0 ? 1 : grid.GetSize()-1;
|
|---|
| 125 | for (Int_t ie=0; ie<nbins; ie++)
|
|---|
| 126 | {
|
|---|
| 127 | switch (ver)
|
|---|
| 128 | {
|
|---|
| 129 | case 0: // Replace Energy Grid by classification
|
|---|
| 130 | {
|
|---|
| 131 | Int_t irows=0;
|
|---|
| 132 | for (Int_t j=0; j<nrows; j++)
|
|---|
| 133 | {
|
|---|
| 134 | const Double_t energy = matrixtrain.GetM()(j,ncols-1);
|
|---|
| 135 | const Bool_t inside = energy>grid[ie] && energy<=grid[ie+1];
|
|---|
| 136 |
|
|---|
| 137 | mat(j, ncols-nobs) = inside ? 1 : 0;
|
|---|
| 138 |
|
|---|
| 139 | if (inside)
|
|---|
| 140 | irows++;
|
|---|
| 141 | }
|
|---|
| 142 | if (irows==0)
|
|---|
| 143 | *fLog << warn << "WARNING - Skipping";
|
|---|
| 144 | else
|
|---|
| 145 | *fLog << inf << "Training RF for";
|
|---|
| 146 |
|
|---|
| 147 | *fLog << " energy bin " << ie << " (" << grid[ie] << ", " << grid[ie+1] << ") " << irows << "/" << nrows << endl;
|
|---|
| 148 |
|
|---|
| 149 | if (irows==0)
|
|---|
| 150 | continue;
|
|---|
| 151 | }
|
|---|
| 152 | break;
|
|---|
| 153 |
|
|---|
| 154 | case 1: // Use Energy as classifier
|
|---|
| 155 | case 2:
|
|---|
| 156 | for (Int_t j=0; j<nrows; j++)
|
|---|
| 157 | mat(j, ncols-nobs) = matrixtrain.GetM()(j,ncols-1);
|
|---|
| 158 | break;
|
|---|
| 159 | }
|
|---|
| 160 |
|
|---|
| 161 | MHMatrix matrix(mat, &rules, "MatrixTrain");
|
|---|
| 162 |
|
|---|
| 163 | MParList plist;
|
|---|
| 164 | MTaskList tlist;
|
|---|
| 165 | plist.AddToList(&tlist);
|
|---|
| 166 | plist.AddToList(&matrix);
|
|---|
| 167 |
|
|---|
| 168 | MRanForest rf;
|
|---|
| 169 | rf.SetNumTrees(fNumTrees);
|
|---|
| 170 | rf.SetNumTry(fNumTry);
|
|---|
| 171 | rf.SetNdSize(fNdSize);
|
|---|
| 172 | rf.SetClassify(ver<2 ? 1 : 0);
|
|---|
| 173 | if (ver==1)
|
|---|
| 174 | rf.SetGrid(grid);
|
|---|
| 175 |
|
|---|
| 176 | plist.AddToList(&rf);
|
|---|
| 177 |
|
|---|
| 178 | MRanForestGrow rfgrow;
|
|---|
| 179 | tlist.AddToList(&rfgrow);
|
|---|
| 180 |
|
|---|
| 181 | MEvtLoop evtloop;
|
|---|
| 182 | evtloop.SetParList(&plist);
|
|---|
| 183 |
|
|---|
| 184 | if (!evtloop.Eventloop())
|
|---|
| 185 | return kFALSE;
|
|---|
| 186 |
|
|---|
| 187 | if (fDebug)
|
|---|
| 188 | gLog.SetNullOutput(kFALSE);
|
|---|
| 189 |
|
|---|
| 190 | if (ver==0)
|
|---|
| 191 | {
|
|---|
| 192 | // Calculate bin center
|
|---|
| 193 | const Double_t E = (TMath::Log10(grid[ie])+TMath::Log10(grid[ie+1]))/2;
|
|---|
| 194 |
|
|---|
| 195 | // save whole forest
|
|---|
| 196 | rf.SetUserVal(E);
|
|---|
| 197 | rf.SetName(Form("%.10f", E));
|
|---|
| 198 | }
|
|---|
| 199 |
|
|---|
| 200 | rf.Write();
|
|---|
| 201 | }
|
|---|
| 202 |
|
|---|
| 203 | // save rules
|
|---|
| 204 | usedrules.Write("rules");
|
|---|
| 205 |
|
|---|
| 206 | return kTRUE;
|
|---|
| 207 | }
|
|---|
| 208 |
|
|---|
| 209 | Int_t MRFEnergyEst::ReadForests(MParList &plist)
|
|---|
| 210 | {
|
|---|
| 211 | TFile fileRF(fFileName, "read");
|
|---|
| 212 | if (!fileRF.IsOpen())
|
|---|
| 213 | {
|
|---|
| 214 | *fLog << err << dbginf << "File containing RFs could not be opened... aborting." << endl;
|
|---|
| 215 | return kFALSE;
|
|---|
| 216 | }
|
|---|
| 217 |
|
|---|
| 218 | fEForests.Delete();
|
|---|
| 219 |
|
|---|
| 220 | TIter Next(fileRF.GetListOfKeys());
|
|---|
| 221 | TObject *o=0;
|
|---|
| 222 | while ((o=Next()))
|
|---|
| 223 | {
|
|---|
| 224 | MRanForest *forest=0;
|
|---|
| 225 | fileRF.GetObject(o->GetName(), forest);
|
|---|
| 226 | if (!forest)
|
|---|
| 227 | continue;
|
|---|
| 228 |
|
|---|
| 229 | forest->SetUserVal(atof(o->GetName()));
|
|---|
| 230 |
|
|---|
| 231 | fEForests.Add(forest);
|
|---|
| 232 | }
|
|---|
| 233 |
|
|---|
| 234 | // Maybe fEForests[0].fRules yould be used instead?
|
|---|
| 235 |
|
|---|
| 236 | if (fData->Read("rules")<=0)
|
|---|
| 237 | {
|
|---|
| 238 | *fLog << err << "ERROR - Reading 'rules' from file " << fFileName << endl;
|
|---|
| 239 | return kFALSE;
|
|---|
| 240 | }
|
|---|
| 241 |
|
|---|
| 242 | return kTRUE;
|
|---|
| 243 | }
|
|---|
| 244 |
|
|---|
| 245 | Int_t MRFEnergyEst::PreProcess(MParList *plist)
|
|---|
| 246 | {
|
|---|
| 247 | fEnergyEst = (MParameterD*)plist->FindCreateObj("MParameterD", "MEnergyEst");
|
|---|
| 248 | if (!fEnergyEst)
|
|---|
| 249 | return kFALSE;
|
|---|
| 250 |
|
|---|
| 251 | fData = (MDataArray*)plist->FindCreateObj("MDataArray");
|
|---|
| 252 | if (!fData)
|
|---|
| 253 | return kFALSE;
|
|---|
| 254 |
|
|---|
| 255 | if (!ReadForests(*plist))
|
|---|
| 256 | {
|
|---|
| 257 | *fLog << err << "Reading RFs failed... aborting." << endl;
|
|---|
| 258 | return kFALSE;
|
|---|
| 259 | }
|
|---|
| 260 |
|
|---|
| 261 | *fLog << inf << "RF read from " << fFileName << endl;
|
|---|
| 262 |
|
|---|
| 263 | if (fTestMatrix)
|
|---|
| 264 | return kTRUE;
|
|---|
| 265 |
|
|---|
| 266 | fData->Print();
|
|---|
| 267 |
|
|---|
| 268 | if (!fData->PreProcess(plist))
|
|---|
| 269 | {
|
|---|
| 270 | *fLog << err << "PreProcessing of the MDataArray failed... aborting." << endl;
|
|---|
| 271 | return kFALSE;
|
|---|
| 272 | }
|
|---|
| 273 |
|
|---|
| 274 | return kTRUE;
|
|---|
| 275 | }
|
|---|
| 276 |
|
|---|
| 277 | #include <TGraph.h>
|
|---|
| 278 | #include <TF1.h>
|
|---|
| 279 | Int_t MRFEnergyEst::Process()
|
|---|
| 280 | {
|
|---|
| 281 | TVector event;
|
|---|
| 282 | if (fTestMatrix)
|
|---|
| 283 | *fTestMatrix >> event;
|
|---|
| 284 | else
|
|---|
| 285 | *fData >> event;
|
|---|
| 286 |
|
|---|
| 287 | // --------------- Single Tree RF -------------------
|
|---|
| 288 | if (fEForests.GetEntries()==1)
|
|---|
| 289 | {
|
|---|
| 290 | MRanForest *rf = (MRanForest*)fEForests[0];
|
|---|
| 291 | fEnergyEst->SetVal(rf->CalcHadroness(event));
|
|---|
| 292 | fEnergyEst->SetReadyToSave();
|
|---|
| 293 |
|
|---|
| 294 | return kTRUE;
|
|---|
| 295 | }
|
|---|
| 296 |
|
|---|
| 297 | // --------------- Multi Tree RF -------------------
|
|---|
| 298 | static TF1 f1("f1", "gaus");
|
|---|
| 299 |
|
|---|
| 300 | Double_t sume = 0;
|
|---|
| 301 | Double_t sumh = 0;
|
|---|
| 302 | Double_t maxh = 0;
|
|---|
| 303 | Double_t maxe = 0;
|
|---|
| 304 |
|
|---|
| 305 | Double_t max = -1e10;
|
|---|
| 306 | Double_t min = 1e10;
|
|---|
| 307 |
|
|---|
| 308 | //TH1C h("", "", fEForests.GetSize(), 0, 1);
|
|---|
| 309 |
|
|---|
| 310 | TIter Next(&fEForests);
|
|---|
| 311 | MRanForest *rf = 0;
|
|---|
| 312 |
|
|---|
| 313 | TGraph g;
|
|---|
| 314 | while ((rf=(MRanForest*)Next()))
|
|---|
| 315 | {
|
|---|
| 316 | const Double_t h = rf->CalcHadroness(event);
|
|---|
| 317 | const Double_t e = rf->GetUserVal();
|
|---|
| 318 |
|
|---|
| 319 | g.SetPoint(g.GetN(), e, h);
|
|---|
| 320 |
|
|---|
| 321 | sume += e*h;
|
|---|
| 322 | sumh += h;
|
|---|
| 323 |
|
|---|
| 324 | if (h>maxh)
|
|---|
| 325 | {
|
|---|
| 326 | maxh = h;
|
|---|
| 327 | maxe = e;
|
|---|
| 328 | }
|
|---|
| 329 | if (e>max)
|
|---|
| 330 | max = e;
|
|---|
| 331 | if (e<min)
|
|---|
| 332 | min = e;
|
|---|
| 333 | }
|
|---|
| 334 |
|
|---|
| 335 | switch (fEstimationMode)
|
|---|
| 336 | {
|
|---|
| 337 | case kMean:
|
|---|
| 338 | fEnergyEst->SetVal(pow(10, sume/sumh));
|
|---|
| 339 | break;
|
|---|
| 340 | case kMaximum:
|
|---|
| 341 | fEnergyEst->SetVal(pow(10, maxe));
|
|---|
| 342 | break;
|
|---|
| 343 | case kFit:
|
|---|
| 344 | f1.SetParameter(0, maxh);
|
|---|
| 345 | f1.SetParameter(1, maxe);
|
|---|
| 346 | f1.SetParameter(2, 0.125);
|
|---|
| 347 | g.Fit(&f1, "Q0N");
|
|---|
| 348 | fEnergyEst->SetVal(pow(10, f1.GetParameter(1)));
|
|---|
| 349 | break;
|
|---|
| 350 | }
|
|---|
| 351 |
|
|---|
| 352 | fEnergyEst->SetReadyToSave();
|
|---|
| 353 |
|
|---|
| 354 | return kTRUE;
|
|---|
| 355 | }
|
|---|
| 356 |
|
|---|
| 357 | // --------------------------------------------------------------------------
|
|---|
| 358 | //
|
|---|
| 359 | //
|
|---|
| 360 | Int_t MRFEnergyEst::ReadEnv(const TEnv &env, TString prefix, Bool_t print)
|
|---|
| 361 | {
|
|---|
| 362 | Bool_t rc = kFALSE;
|
|---|
| 363 | if (IsEnvDefined(env, prefix, "FileName", print))
|
|---|
| 364 | {
|
|---|
| 365 | rc = kTRUE;
|
|---|
| 366 | SetFileName(GetEnvValue(env, prefix, "FileName", fFileName));
|
|---|
| 367 | }
|
|---|
| 368 | if (IsEnvDefined(env, prefix, "Debug", print))
|
|---|
| 369 | {
|
|---|
| 370 | rc = kTRUE;
|
|---|
| 371 | SetDebug(GetEnvValue(env, prefix, "Debug", fDebug));
|
|---|
| 372 | }
|
|---|
| 373 | if (IsEnvDefined(env, prefix, "EstimationMode", print))
|
|---|
| 374 | {
|
|---|
| 375 | TString txt = GetEnvValue(env, prefix, "EstimationMode", "");
|
|---|
| 376 | txt = txt.Strip(TString::kBoth);
|
|---|
| 377 | txt.ToLower();
|
|---|
| 378 | if (txt==(TString)"mean")
|
|---|
| 379 | fEstimationMode = kMean;
|
|---|
| 380 | if (txt==(TString)"maximum")
|
|---|
| 381 | fEstimationMode = kMaximum;
|
|---|
| 382 | if (txt==(TString)"fit")
|
|---|
| 383 | fEstimationMode = kFit;
|
|---|
| 384 | rc = kTRUE;
|
|---|
| 385 | }
|
|---|
| 386 | return rc;
|
|---|
| 387 | }
|
|---|